00001 
00002 
00003 
00004 
00005 
00006 
00007 
00008 
00009 
00010 
00011 
00012 
00013       subroutine prism_enddef (ierror)
00014 
00015 
00016 
00017 
00018       use PRISM, dummy_interface => prism_enddef
00019       use psmile_user_data, only : psmile_apply_user_data, psmile_merge_fields
00020       use PSMILe
00021       use PSMILe_SMIOC, only : sga_smioc_comp, transient
00022       use psmile_timer, only : psmile_timer_stop
00023 
00024       implicit none
00025 
00026 
00027 
00028       integer, Intent (Out)               :: ierror
00029 
00030 
00031 
00032 
00033 
00034 
00035 
00036 
00037 
00038 
00039 
00040 
00041 
00042 
00043 
00044       integer, parameter           :: nd_acomps = 3
00045       integer                      :: a_comps (nd_acomps, noComponents)
00046       Integer                      :: n_active
00047       Integer                      :: my_icomp0_coll_comps
00048       Integer                      :: icomp, lastc
00049       Integer                      :: use_Grid(Number_of_Grids_allocated)
00050 
00051       Type (Transient), pointer    :: sga_smioc_transi(:)
00052       Logical                      :: found
00053 
00054       Integer                      :: i, j, n
00055 
00056       logical, save                :: first = .true.
00057 
00058 
00059 
00060       integer                      :: nbr_in, nbr_out, il_side, il_i, il_o
00061       integer                      :: il_interp_meth
00062       integer                      :: il_userdef_id, il_ass_var_id
00063 
00064 
00065 
00066       integer                      :: ninter, nmyint, nnull, max_num_intersect
00067       integer, allocatable         :: num_intersect_per_grid(:), num_dummy_intersect_per_grid(:)
00068 
00069 
00070 
00071       integer                      :: lastag
00072       integer                      :: grid_id
00073       logical                      :: get_halo
00074       integer                      :: color
00075       integer                      :: key
00076 
00077 
00078 
00079       Integer, parameter           :: nerrp = 3
00080       Integer                      :: ierrp (nerrp)
00081 
00082 #ifdef TIMING
00083 
00084 
00085       Double Precision            :: tic, toc
00086 #endif
00087 
00088 
00089 
00090 
00091 
00092 
00093 
00094 
00095 
00096 
00097 
00098 
00099 
00100 
00101 
00102 
00103 
00104 
00105 
00106 
00107 
00108   Character(len=len_cvs_string), save :: mycvs = 
00109       '$Id: prism_enddef.F90 3248 2011-06-23 13:03:19Z coquart $'
00110 
00111 
00112 
00113 #ifdef VERBOSE
00114       print 9990, trim(ch_id)
00115       call psmile_flushstd
00116 #endif /* VERBOSE */
00117 
00118       if ( .not. first ) then
00119          ierror = -1
00120 #ifdef VERBOSE
00121          print *, trim(ch_id), ': prism_enddef: It is not allowed to call prism_enddef twice.'
00122          print *, trim(ch_id), ': prism_enddef: eof ierror =', ierror
00123 #endif /* VERBOSE */
00124          return
00125       endif
00126 
00127 #ifdef TIMING
00128       call MPI_Barrier ( comm_psmile, ierror )
00129       tic=MPI_Wtime()
00130 #endif
00131 
00132 
00133 
00134       ierror = 0
00135       lastag = PSMILe_Enddef_Tag
00136       first  = .false.
00137       use_Grid = PSMILe_Undef
00138       get_halo = .false.
00139 
00140 
00141 
00142       call psmile_apply_user_data(ierror)
00143 
00144 
00145 
00146 
00147 #ifdef DEBUG
00148             print *, trim(ch_id),': prism_enddef: Fields allocated = ',Number_of_Fields_allocated
00149 #endif
00150       do i = 1, Number_of_Fields_allocated
00151 
00152          if ( Fields(i)%status  == PSMILe_status_defined .and. &
00153               Fields(i)%used_for_coupling ) then
00154 
00155 
00156 
00157 
00158 
00159             nbr_in  = 0
00160             nbr_out = 0
00161 
00162             nbr_in = Fields(i)%Taskin%nbr_inchannels
00163             if ( Associated(Fields(i)%Taskout) ) then
00164                nbr_out = size (Fields(i)%Taskout)
00165             endif
00166 #ifdef DEBUG
00167             print *, trim(ch_id),': prism_enddef: index Field = ',i
00168             print *, trim(ch_id),': prism_enddef: nbr_in nbr_out = ',nbr_in,nbr_out
00169 #endif /* DEBUG */
00170 
00171 
00172 
00173 
00174 
00175          do il_i = 1, nbr_in
00176             il_side = 1
00177 
00178             Fields(i)%Taskin%In_channel(il_i)%assoc_var_id  = PSMILe_undef
00179             Fields(i)%Taskin%In_channel(il_i)%userdef_id    = PSMILe_undef
00180             Grids(Methods(Fields(i)%method_id)%grid_id)%assoc_grid_id = PSMILe_undef
00181 
00182             il_interp_meth = Fields(i)%Taskin%In_channel(il_i)%interp%interp_meth(1)
00183             if ( il_interp_meth == PSMILe_user3D ) then
00184 
00185 
00186                call psmile_set_userdef(i, il_side, il_i, ierror)
00187 #ifdef DEBUG
00188       print *, trim(ch_id), ': prism_enddef: index Field = ',i
00189       print *, trim(ch_id), ': prism_enddef: index channel in ',il_i
00190       print *, trim(ch_id), ': prism_enddef: assoc_var_id = ',  &
00191                                Fields(i)%Taskin%In_channel(il_i)%assoc_var_id
00192       print *, trim(ch_id), ': prism_enddef: userdef_id = ',    &
00193                                Fields(i)%Taskin%In_channel(il_i)%userdef_id
00194       print *, trim(ch_id), ': prism_enddef: assoc_grid_id = ', &
00195                                Grids(Methods(Fields(i)%method_id)%grid_id)%assoc_grid_id
00196 #endif /* DEBUG */
00197             endif
00198          enddo        
00199 
00200 
00201          do il_o = 1, nbr_out
00202             il_side = 0
00203 
00204             Fields(i)%Taskout(il_o)%assoc_var_id  = PSMILe_undef
00205             Fields(i)%Taskout(il_o)%userdef_id    = PSMILe_undef
00206             Grids(Methods(Fields(i)%method_id)%grid_id)%assoc_grid_id = PSMILe_undef
00207 
00208             il_interp_meth = Fields(i)%Taskout(il_o)%interp%interp_meth(1)
00209             if ( il_interp_meth == PSMILe_user3D ) then
00210 
00211 
00212                call psmile_set_userdef(i, il_side, il_o, ierror)
00213 #ifdef DEBUG
00214       print *, trim(ch_id), ': prism_enddef: index Field = ',i
00215       print *, trim(ch_id), ': prism_enddef: index channel out ',il_o
00216       print *, trim(ch_id), ': prism_enddef: assoc_var_id = ',  &
00217                                Fields(i)%Taskout(il_o)%assoc_var_id
00218       print *, trim(ch_id), ': prism_enddef: userdef_id = ',    &
00219                                Fields(i)%Taskout(il_o)%userdef_id
00220       print *, trim(ch_id), ': prism_enddef: assoc_grid_id = ', &
00221                                Grids(Methods(Fields(i)%method_id)%grid_id)%assoc_grid_id
00222 #endif /* DEBUG */
00223             endif
00224          enddo                
00225 
00226          endif             
00227 
00228       end do            
00229 
00230 
00231 
00232       call psmile_control_grids (ierror)
00233       if (ierror /= 0) return
00234 
00235 
00236 
00237       do n = 1, Number_of_Comps_allocated
00238          if ( Comps(n)%status == PSMILe_status_defined ) then
00239             sga_smioc_transi => sga_smioc_comp(n)%sga_smioc_transi
00240 
00241             do i = 1, size(sga_smioc_transi)
00242 
00243                found = .false.
00244                do j = 1, Number_of_Fields_allocated
00245                   if ( Fields(j)%status == PSMILe_status_defined ) then
00246                      if ( trim(sga_smioc_transi(i)%cg_local_name) == trim(adjustl(Fields(j)%local_name)) ) &
00247                         found = .true.
00248                   endif
00249                enddo
00250                if ( .not. found ) then
00251                   ierror =  PRISM_Warn_NoDefVar
00252                   ierrp (1) = j
00253                   print *, trim(ch_id), ': transient with name ', &
00254                        trim(sga_smioc_transi(i)%cg_local_name), ' failed.'
00255                    call psmile_warning ( ierror,  trim(sga_smioc_transi(i)%cg_local_name), &
00256                         ierrp(1), 1, __FILE__, __LINE__ )
00257                endif
00258 
00259             enddo
00260 
00261          endif
00262       enddo
00263 
00264 
00265 
00266 
00267       call psmile_get_act_comps (a_comps, nd_acomps, n_act_comp, ierror)
00268       if (ierror /= 0) return
00269 
00270 
00271 
00272       Allocate (comp_infos(1:n_act_comp), STAT = ierror)
00273       if ( ierror > 0 ) then
00274          ierrp (1) = ierror
00275          ierrp (2) = n_act_comp
00276          call psmile_error ( PRISM_Error_Alloc, 'comp_infos', &
00277                              ierrp, 2, __FILE__, __LINE__ )
00278          return
00279       endif
00280 
00281 
00282 
00283 
00284 
00285 
00286 
00287       do grid_id = 1, Number_of_Grids_allocated
00288       if ( Grids(grid_id)%grid_type == PRISM_Gaussreduced_regvrt ) then
00289          call psmile_gauss_setup ( grid_id, ierror )
00290          if (ierror /= 0) return
00291       endif
00292       enddo
00293 
00294 
00295 
00296          do i = 1, n_act_comp
00297          call psmile_enddef_comp (a_comps(3,i), a_comps (1,i), &
00298                                   a_comps(2,i), comp_infos(i), ierror)
00299          if (ierror /= 0) return
00300          end do
00301 
00302 
00303 
00304 
00305 
00306 
00307 
00308 
00309 
00310 
00311 
00312       call psmile_enddef_appl (lastag, my_icomp0_coll_comps, &
00313                                n_active, ierror)
00314       if (ierror /= 0) return
00315 
00316 #ifdef TIMING
00317       call MPI_Barrier ( comm_psmile, ierror )
00318       toc=MPI_Wtime()
00319       print *, trim(ch_id), ': prism_enddef: Time for initial procedures ', toc-tic
00320 #endif
00321 
00322 #ifdef __PSMILE_WITH_IO
00323 
00324 
00325 
00326 
00327 #ifdef PRISM_ASSERTION
00328       IF (PRISM_noCompsPerProc /= n_act_comp) THEN
00329           CALL psmile_assert ( __FILE__, __LINE__, &
00330                 'Call to psmile_io_init_pelist supposes that PRISM_noCompsPerProc = n_act_comp')
00331       ENDIF
00332 #endif
00333       do i = 1, PRISM_noCompsPerProc
00334          call psmile_io_init_pelist (i, comp_infos(i),  ierror)
00335          if (ierror /= 0) return         
00336       enddo
00337 
00338       call psmile_enddef_metadata(ierror)
00339 
00340       if ( ierror .ne. 0 ) then
00341         ierrp(1)=ierror
00342         call psmile_error(PRISM_Error_IO_Meta, &
00343                          'psmile_enddef_metadata', &
00344                               ierrp, 1, __FILE__, __LINE__ )
00345       endif
00346 
00347       if ( Appl%stand_alone ) then
00348 #ifdef VERBOSE
00349          print *, trim(ch_id), ': prism_enddef: Stand alone application. We return'
00350          print *, trim(ch_id), ': prism_enddef: eof ierror =', ierror
00351 #endif /* VERBOSE */
00352          return
00353       endif
00354 
00355 
00356 #endif
00357 
00358 #ifdef TIMING
00359       call MPI_Barrier ( comm_psmile, ierror )
00360       tic=MPI_Wtime()
00361 #endif
00362 
00363 
00364 
00365 
00366       ninter = 0
00367       nmyint = 0
00368       nnull  = 0
00369 
00370       lastc = 0
00371 
00372 #define CONS_REMAP_DEADLOCK_WORKAROUND
00373 #ifdef CONS_REMAP_DEADLOCK_WORKAROUND
00374 
00375       max_num_intersect = 1
00376       do i = 1,  Number_of_coll_comps
00377 
00378          
00379          
00380          
00381          max_num_intersect = max_num_intersect * &
00382             SUM(all_comp_infos(i)%Number_of_grids_vector(1:all_comp_infos(i)%size))
00383       enddo
00384 
00385       call psmile_flushstd
00386       allocate (paction%intersect_ranks(max_num_intersect), stat = ierror)
00387       if ( ierror > 0 ) then
00388          ierrp (1) = ierror
00389          ierrp (2) = max_num_intersect
00390          call psmile_error ( PRISM_Error_Alloc, 'intersect_ranks', &
00391                              ierrp, 2, __FILE__, __LINE__ )
00392          return
00393       endif
00394       paction%intersect_ranks = -1
00395 #endif
00396 
00397       allocate (num_intersect_per_grid(Number_of_Grids_allocated), &
00398                 num_dummy_intersect_per_grid(Number_of_Grids_allocated), stat = ierror)
00399       if ( ierror > 0 ) then
00400          ierrp (1) = ierror
00401          ierrp (2) = Number_of_Grids_allocated
00402          call psmile_error ( PRISM_Error_Alloc, 'num_intersect_per_grid, num_dummy_intersect_per_grid', &
00403                              ierrp, 2, __FILE__, __LINE__ )
00404          return
00405       endif
00406 
00407       num_intersect_per_grid(:)        = 0
00408       num_dummy_intersect_per_grid (:) = 0
00409 
00410       do i = 1, n_act_comp
00411 
00412 
00413 
00414         do icomp = lastc+1, n_active
00415            if (all_comp_infos(my_icomp0_coll_comps+icomp)%global_comp_id ==   &
00416                                             comp_infos(i)%global_comp_id) exit
00417         end do
00418 
00419         if (icomp > n_active) then
00420            ierror = PRISM_Error_Internal
00421            ierrp (1) = comp_infos(i)%global_comp_id
00422            ierrp (2) = lastc
00423            ierrp (3) = n_active
00424 
00425            call psmile_error ( ierror, 'Global Comp Id not found in all_comp_infos', &
00426                                ierrp, 3, __FILE__, __LINE__ )
00427         endif
00428 
00429         lastc = icomp
00430 
00431 
00432 
00433         call psmile_find_intersect (comp_infos(i), my_icomp0_coll_comps+icomp, &
00434                                     num_intersect_per_grid, &
00435                                     num_dummy_intersect_per_grid, &
00436                                     ninter, nmyint, nnull, lastag, ierror)
00437         if (ierror /= 0) return
00438       end do
00439 
00440 
00441 
00442 
00443 
00444       color = 0
00445       key   = 0
00446 
00447       do n = 1, Number_of_Comps_allocated
00448 
00449         if ( Comps(n)%status /= PSMILe_status_free ) then
00450 
00451         do i = 1, Number_of_Fields_allocated
00452            if ( Fields(i)%smioc_loc /= PSMILe_Undef .and. Fields(i)%comp_id == n ) then
00453               if ( associated(Fields(i)%Taskout) ) then
00454                  grid_id = Methods(Fields(i)%method_id)%grid_id
00455                  if ( Associated (Grids(grid_id)%partition) ) then
00456                     if ( ( Grids(grid_id)%grid_type == PRISM_Irrlonlat_regvrt &
00457 #ifdef TODO
00458 
00459                       .or. Grids(grid_id)%grid_type == PRISM_Irrlonlatvrt
00460 #endif
00461                          )  .and. &
00462                          use_Grid(grid_id) == PSMILe_Undef ) then
00463                          use_Grid(grid_id) = Grids(grid_id)%global_grid_id
00464                          get_halo = .true.
00465                          color    = 1
00466                     endif
00467                  endif
00468               endif
00469            endif
00470         enddo
00471 
00472         call MPI_Comm_split (Comps(n)%comm, color, key, Comps(n)%comm_halo, ierror)
00473 
00474         if ( get_halo ) then
00475 
00476            call psmile_get_halo_indices ( n, use_Grid, ierror )
00477 
00478            call psmile_get_halo_points ( n, ierror )
00479 
00480         endif
00481 
00482         endif
00483 
00484       enddo 
00485 
00486 
00487 
00488 
00489 
00490 
00491       call psmile_get_intersect (ninter, nmyint, nnull, &
00492                                  num_intersect_per_grid, &
00493                                  num_dummy_intersect_per_grid, &
00494                                  lastag, ierror)
00495       if (ierror /= 0) return
00496 
00497 #ifdef TIMING
00498       call MPI_Barrier ( comm_psmile, ierror )
00499       toc=MPI_Wtime()
00500       print *, trim(ch_id), ': prism_enddef: Time for search', toc-tic
00501 #endif
00502 
00503 
00504 
00505 
00506 
00507       call psmile_mg_clean (ierror)
00508       if (ierror /= 0) return
00509 
00510 
00511 
00512 
00513 
00514 
00515 
00516 
00517 
00518 
00519 #define CONS_REMAP_DEADLOCK_WORKAROUND
00520 #ifdef CONS_REMAP_DEADLOCK_WORKAROUND
00521       Deallocate (paction%intersect_ranks)
00522 #endif
00523       Deallocate (num_intersect_per_grid,num_dummy_intersect_per_grid)
00524       Deallocate (comp_infos)
00525 
00526       if (Number_of_coll_comps > 0) then
00527          do i = 1, Number_of_coll_comps
00528             Deallocate (all_comp_infos(i)%Number_of_Grids_vector)
00529             Deallocate (all_comp_infos(i)%all_extent_infos)
00530             Deallocate (all_comp_infos(i)%psmile_ranks)
00531          end do
00532 
00533          Deallocate (all_comp_infos)
00534          Number_of_coll_comps = 0
00535       endif
00536 
00537 
00538 
00539 
00540 
00541 
00542       call psmile_get_restart ( ierror )
00543 
00544       if ( ierror .ne. 0 ) then
00545         ierrp(1)=ierror
00546         call psmile_error(999, &
00547                          'psmile_get_restart', &
00548                               ierrp, 1, __FILE__, __LINE__ )
00549       endif
00550 
00551 #ifdef DEBUG
00552 
00553 
00554 
00555       do i = 1, Number_of_Grids_allocated
00556         call psmile_print_grid_info (i)
00557       enddo
00558 
00559       do i = 1, Number_of_Fields_allocated
00560         call psmile_print_field_info (i)
00561       enddo
00562 
00563       do i = 1, Number_of_Methods_allocated
00564         call psmile_print_method_info (i)
00565       enddo
00566 #endif
00567 
00568 
00569 
00570 
00571 
00572       call psmile_merge_fields(ierror)
00573 
00574 
00575 
00576 #ifdef VERBOSE
00577       print 9980, trim(ch_id), ierror
00578       call psmile_flushstd
00579 #endif /* VERBOSE */
00580 
00581 
00582 
00583 #ifdef PROFILE
00584       call psmile_timer_stop(2)
00585 #endif
00586 
00587 
00588 
00589 #ifdef VERBOSE
00590 9990 format (1x, a, ': prism_enddef: start')
00591 9980 format (1x, a, ': prism_enddef: eof ierror =', i3)
00592 #endif /* VERBOSE */
00593 
00594       end subroutine prism_enddef