55 integer,
intent(in) :: dir
56 double precision,
intent(in) :: xslice
58 double precision,
dimension(0:nw+nwauxio) :: normconv
59 integer :: Njgrid, jgrid
60 integer,
dimension(ndim-1) :: ixsubGlo, ixsubGhi
61 integer,
dimension(ndim-1) :: ixsubMlo, ixsubMhi
62 integer :: size_subblock_io, nx^D, slice_fh, nwexpand
63 integer :: type_subblock_io, type_subblockC_io, type_subblock_x_io, type_subblockC_x_io
64 integer,
dimension(ndim) :: sizes, subsizes, start
67 nx^d=ixmhi^d-ixmlo^d+1;
83 if(xslice<xprobmin^d.or.xslice>xprobmax^d) &
84 call mpistop(
"slice out of bounds")
95 sizes(1) = ixghi2; sizes(2) = ixghi3;
96 ixsubglo(1) = ixglo2; ixsubglo(2) = ixglo3;
97 ixsubghi(1) = ixghi2; ixsubghi(2) = ixghi3;
98 subsizes(1)=nx2;subsizes(2)=nx3;
99 start(1)=ixmlo2-1;start(2)=ixmlo3-1;
102 sizes(1) = ixghi1; sizes(2) = ixghi3;
103 ixsubglo(1) = ixglo1; ixsubglo(2) = ixglo3;
104 ixsubghi(1) = ixghi1; ixsubghi(2) = ixghi3;
105 subsizes(1)=nx1;subsizes(2)=nx3;
106 start(1)=ixmlo1-1;start(2)=ixmlo3-1;
109 ixsubglo(1) = ixglo1; ixsubglo(2) = ixglo2;
110 ixsubghi(1) = ixghi1; ixsubghi(2) = ixghi2;
111 sizes(1) = ixghi1; sizes(2) = ixghi2;
112 subsizes(1)=nx1;subsizes(2)=nx2;
113 start(1)=ixmlo1-1;start(2)=ixmlo2-1;
116 call mpistop(
"slice direction not clear in put_slice")
122 ixsubglo(1) = ixglo2; ixsubghi(1) = ixghi2;
128 ixsubglo(1) = ixglo1; ixsubghi(1) = ixghi1;
134 call mpistop(
"slice direction not clear in put_slice")
142 {^de&ixsubmlo(^de-1) = ixsubglo(^de-1)+
nghostcells;}
143 {^de&ixsubmhi(^de-1) = ixsubghi(^de-1)-
nghostcells;}
146 sizes(
ndim)=nw+nwexpand
147 subsizes(
ndim)=nw+nwexpand
151 call mpi_type_create_subarray(
ndim,sizes,subsizes,start, &
152 mpi_order_fortran,mpi_double_precision, &
154 call mpi_type_commit(type_subblock_io,
ierrmpi)
159 call mpi_type_create_subarray(
ndim,sizes,subsizes,start, &
160 mpi_order_fortran,mpi_double_precision, &
162 call mpi_type_commit(type_subblock_x_io,
ierrmpi)
166 subsizes(1:
ndim-1) = subsizes(1:
ndim-1) + 1
167 start(1:
ndim-1) = start(1:
ndim-1) - 1
168 sizes(
ndim)=nw+nwexpand
169 subsizes(
ndim)=nw+nwexpand
172 call mpi_type_create_subarray(
ndim,sizes,subsizes,start, &
173 mpi_order_fortran,mpi_double_precision, &
175 call mpi_type_commit(type_subblockc_io,
ierrmpi)
180 call mpi_type_create_subarray(
ndim,sizes,subsizes,start, &
181 mpi_order_fortran,mpi_double_precision, &
183 call mpi_type_commit(type_subblockc_x_io,
ierrmpi)
198 case (
'vtu',
'vtuCC')
207 call mpi_type_free(type_subblock_io,
ierrmpi)
208 call mpi_type_free(type_subblock_x_io,
ierrmpi)
209 call mpi_type_free(type_subblockc_io,
ierrmpi)
210 call mpi_type_free(type_subblockc_x_io,
ierrmpi)
218 integer :: status(MPI_STATUS_SIZE), ipe
220 character(len=1024) :: filename, xlabel
221 character(len=79) :: xxlabel
222 character(len=name_len) :: wnamei(1:nw+nwauxio),xandwnamei(1:ndim+nw+nwauxio)
223 character(len=1024) :: outfilehead
227 inquire(slice_fh,opened=fileopen)
228 if(.not.fileopen)
then
230 write(xlabel,
"(D9.2)")xslice
233 write(xxlabel(1:1),
"(a)")
"+"
235 write(filename,
"(a,i1.1,a,i4.4,a)") trim(
base_filename)//
'_d',dir,
'_x'//trim(xxlabel)//
'_n',
slicenext,
'.vtu'
236 open(slice_fh,file=filename,status=
'unknown',form=
'formatted')
241 write(slice_fh,
'(a)')
'<?xml version="1.0"?>'
242 write(slice_fh,
'(a)',advance=
'no')
'<VTKFile type="UnstructuredGrid"'
244 write(slice_fh,
'(a)')
' version="0.1" byte_order="LittleEndian">'
246 write(slice_fh,
'(a)')
' version="0.1" byte_order="BigEndian">'
248 write(slice_fh,
'(a)')
' <UnstructuredGrid>'
249 write(slice_fh,
'(a)')
'<FieldData>'
250 write(slice_fh,
'(2a)')
'<DataArray type="Float32" Name="TIME" ',&
251 'NumberOfTuples="1" format="ascii">'
253 write(slice_fh,
'(a)')
'</DataArray>'
254 write(slice_fh,
'(a)')
'</FieldData>'
258 call write_slice_vtk(jgrid,slice_fh,wnamei)
272 if (ipe ==
mype )
then
273 call mpi_send(ps_sub(jgrid)%x,1,type_subblock_x_io,0,itag,
icomm,
ierrmpi)
274 call mpi_send(ps_sub(jgrid)%w,1,type_subblock_io,0,itag+1,
icomm,
ierrmpi)
275 call mpi_send(ps_sub(jgrid)%xC,1,type_subblockc_x_io,0,itag+2,
icomm,
ierrmpi)
276 call mpi_send(ps_sub(jgrid)%wC,1,type_subblockc_io,0,itag+3,
icomm,
ierrmpi)
277 call mpi_send(normconv,nw+nwauxio+1,mpi_double_precision,0,itag+4,
icomm,
ierrmpi)
280 call mpi_recv(ps_sub(njgrid+1)%x,1,type_subblock_x_io,ipe,itag,
icomm,status,
ierrmpi)
281 call mpi_recv(ps_sub(njgrid+1)%w,1,type_subblock_io,ipe,itag+1,
icomm,status,
ierrmpi)
282 call mpi_recv(ps_sub(njgrid+1)%xC,1,type_subblockc_x_io,ipe,itag+2,
icomm,status,
ierrmpi)
283 call mpi_recv(ps_sub(njgrid+1)%wC,1,type_subblockc_io,ipe,itag+3,
icomm,status,
ierrmpi)
284 call mpi_recv(normconv,nw+nwauxio+1,mpi_double_precision,ipe,&
286 call write_slice_vtk(njgrid+1,slice_fh,wnamei)
294 write(slice_fh,
'(a)')
'</UnstructuredGrid>'
295 write(slice_fh,
'(a)')
'</VTKFile>'
303 subroutine write_slice_vtk(jgrid,slice_fh,wnamei)
306 integer,
intent(in) :: jgrid, slice_fh
307 character(len=name_len),
intent(in) :: wnamei(1:nw+nwauxio)
311 double precision :: x_VTK(1:3)
312 double precision,
parameter :: minvalue = 1.0d-99, maxvalue = 1.0d+99
313 integer :: ixC^L, ixCC^L, nc, np, iw
314 integer :: nx^DM, nxC^DM, icell, ix^DM
317 {^dm&ixccmin^dm = ixsubmlo(^dm);}
318 {^dm&ixccmax^dm = ixsubmhi(^dm);}
319 {^dm&ixcmin^dm = ixsubmlo(^dm)-1;}
320 {^dm&ixcmax^dm = ixsubmhi(^dm);}
322 nx^dm=ixccmax^dm-ixccmin^dm+1;
328 write(slice_fh,
'(a,i7,a,i7,a)') &
329 '<Piece NumberOfPoints="',np,
'" NumberOfCells="',nc,
'">'
337 write(slice_fh,
'(a)')
'<CellData>'
342 write(slice_fh,
'(a,a,a)')&
343 '<DataArray type="Float64" Name="',trim(wnamei(iw)),
'" format="ascii">'
344 write(slice_fh,
'(200(1pe14.6))') {^dm&(|}
roundoff_minmax(ps_sub(jgrid)%w(ix^dm,iw)*normconv(iw),minvalue,maxvalue),{ix^dm=ixccmin^dm,ixccmax^dm)}
345 write(slice_fh,
'(a)')
'</DataArray>'
347 write(slice_fh,
'(a)')
'</CellData>'
351 write(slice_fh,
'(a)')
'<PointData>'
356 write(slice_fh,
'(a,a,a)')&
357 '<DataArray type="Float64" Name="',trim(wnamei(iw)),
'" format="ascii">'
358 write(slice_fh,
'(200(1pe14.6))') {^dm&(|}
roundoff_minmax(ps_sub(jgrid)%wC(ix^dm,iw)*normconv(iw),minvalue,maxvalue),{ix^dm=ixcmin^dm,ixcmax^dm)}
359 write(slice_fh,
'(a)')
'</DataArray>'
361 write(slice_fh,
'(a)')
'</PointData>'
372 write(slice_fh,
'(a)')
'<Points>'
373 write(slice_fh,
'(a)')
'<DataArray type="Float32" NumberOfComponents="3" format="ascii">'
375 {
do ix^dmb=ixcmin^dmb,ixcmax^dmb \}
377 x_vtk(1:
ndim)=ps_sub(jgrid)%xC(ix^dm,1:
ndim)*normconv(0);
378 write(slice_fh,
'(3(1pe14.6))') x_vtk
380 write(slice_fh,
'(a)')
'</DataArray>'
381 write(slice_fh,
'(a)')
'</Points>'
389 write(slice_fh,
'(a)')
'<Cells>'
392 write(slice_fh,
'(a)')
'<DataArray type="Int32" Name="connectivity" format="ascii">'
394 {^dm&
do ix^dmb=1,nx^dmb\}
396 write(slice_fh,
'(2(i7,1x))')ix1-1,ix1
398 write(slice_fh,
'(4(i7,1x))')(ix2-1)*nxc1+ix1-1, &
399 (ix2-1)*nxc1+ix1,ix2*nxc1+ix1-1,ix2*nxc1+ix1
402 write(slice_fh,
'(a)')
'</DataArray>'
405 write(slice_fh,
'(a)')
'<DataArray type="Int32" Name="offsets" format="ascii">'
407 write(slice_fh,
'(i7)') icell*(2**(^nd-1))
409 write(slice_fh,
'(a)')
'</DataArray>'
412 write(slice_fh,
'(a)')
'<DataArray type="Int32" Name="types" format="ascii">'
414 {^iftwod vtk_type=3 \}
415 {^ifthreed vtk_type=8 \}
417 write(slice_fh,
'(i2)') vtk_type
419 write(slice_fh,
'(a)')
'</DataArray>'
421 write(slice_fh,
'(a)')
'</Cells>'
425 write(slice_fh,
'(a)')
'</Piece>'
428 end subroutine write_slice_vtk
430 subroutine put_slice_csv
433 integer :: iw, ipe, itag
434 integer :: status(MPI_STATUS_SIZE)
436 character(len=1024) :: filename, xlabel
437 character(len=79) :: xxlabel
438 character(len=name_len) :: wnamei(1:nw+nwauxio),xandwnamei(1:ndim+nw+nwauxio)
439 character(len=1024) :: outfilehead
440 character(len=1024) :: line
443 inquire(slice_fh,opened=fileopen)
444 if(.not.fileopen)
then
446 write(xlabel,
"(D9.2)")xslice
449 write(xxlabel(1:1),
"(a)")
"+"
451 write(filename,
"(a,i1.1,a,i4.4,a)") trim(base_filename)//
'_d',dir,
'_x'//trim(xxlabel)//
'_n',
slicenext,
'.csv'
452 open(slice_fh,file=filename,status=
'unknown',form=
'formatted')
457 do iw=1,ndim+nw+nwauxio-1
458 line = trim(line)//trim(xandwnamei(iw))//
', '
460 line = trim(line)//trim(xandwnamei(ndim+nw+nwauxio))
461 write(slice_fh,
'(a)')trim(line)
467 call put_slice_line(jgrid,slice_fh)
474 do jgrid=1,morton_sub_stop(ipe)-morton_sub_start(ipe)+1
475 itag=morton_sub_start(ipe)+jgrid-1
476 if (ipe == mype )
then
477 call mpi_send(ps_sub(jgrid)%x,1,type_subblock_x_io,0,itag,icomm,ierrmpi)
478 call mpi_send(ps_sub(jgrid)%w,1,type_subblock_io,0,itag,icomm,ierrmpi)
479 call mpi_send(normconv,nw+nwauxio+1,mpi_double_precision,0,itag,icomm,ierrmpi)
482 call mpi_recv(ps_sub(njgrid+1)%x,1,type_subblock_x_io,ipe,itag,icomm,status,ierrmpi)
483 call mpi_recv(ps_sub(njgrid+1)%w,1,type_subblock_io,ipe,itag,icomm,status,ierrmpi)
484 call mpi_recv(normconv,nw+nwauxio+1,mpi_double_precision,ipe,&
485 itag,icomm,status,ierrmpi)
486 call put_slice_line(njgrid+1,slice_fh)
497 end subroutine put_slice_csv
499 subroutine put_slice_line(jout,file_handle)
500 integer,
intent(in) :: jout, file_handle
502 double precision,
parameter :: minvalue = 1.0d-99, maxvalue = 1.0d+99
503 integer :: ix^D,idir,iw
504 character(len=1024) ::line, data
507 do ix2=ixsubmlo(2),ixsubmhi(2)
508 do ix1=ixsubmlo(1),ixsubmhi(1)
511 do ix1=ixsubmlo(1),ixsubmhi(1)
517 write(
data,
"(es14.6)")
roundoff_minmax(ps_sub(jout)%x(ix1,ix2,idir),minvalue,maxvalue)
520 write(
data,
"(es14.6)")
roundoff_minmax(ps_sub(jout)%x(ix1,idir),minvalue,maxvalue)
523 write(
data,
"(es14.6)")
roundoff_minmax(ps_sub(jout)%x(idir),minvalue,maxvalue)
526 line = trim(line)//trim(data)//
', '
528 do iw = 1,nw+nwauxio-1
530 write(
data,
"(es14.6)")
roundoff_minmax(ps_sub(jout)%w(ix1,ix2,iw)*normconv(iw),minvalue,maxvalue)
533 write(
data,
"(es14.6)")
roundoff_minmax(ps_sub(jout)%w(ix1,iw)*normconv(iw),minvalue,maxvalue)
536 write(
data,
"(es14.6)")
roundoff_minmax(ps_sub(jout)%w(iw)*normconv(iw),minvalue,maxvalue)
538 line = trim(line)//trim(data)//
', '
541 write(
data,
"(es14.6)")
roundoff_minmax(ps_sub(jout)%w(ix1,ix2,nw+nwauxio)*normconv(nw+nwauxio),minvalue,maxvalue)
544 write(
data,
"(es14.6)")
roundoff_minmax(ps_sub(jout)%w(ix1,nw+nwauxio)*normconv(nw+nwauxio),minvalue,maxvalue)
546 line = trim(line)//trim(data)
547 write(file_handle,
'(a)')trim(line)
556 end subroutine put_slice_line
558 subroutine put_slice_dat
560 integer,
dimension(max_blocks) :: iorequest
561 integer,
dimension(MPI_STATUS_SIZE,max_blocks) :: iostatus
562 integer(kind=MPI_OFFSET_KIND) :: offset
564 integer :: amode, status(MPI_STATUS_SIZE), iwrite
565 character(len=1024) :: filename, xlabel
566 character(len=79) :: xxlabel
568 nsubleafs=morton_sub_stop(npe-1)
570 write(xlabel,
"(D9.2)")xslice
573 write(xxlabel(1:1),
"(a)")
"+"
575 write(filename,
"(a,i1.1,a,i4.4,a)") trim(base_filename)//
'_d',dir,
'_x'//trim(xxlabel)//
'_n',
slicenext,
'.dat'
578 open(unit=slice_fh,file=filename,status=
'replace')
582 amode=ior(mpi_mode_create,mpi_mode_wronly)
583 call mpi_file_open(icomm,filename,amode,mpi_info_null,slice_fh,ierrmpi)
584 iorequest=mpi_request_null
589 offset=int(size_subblock_io,kind=mpi_offset_kind) &
590 *int(morton_sub_start(mype)+jgrid-2,kind=mpi_offset_kind)
591 call mpi_file_iwrite_at(slice_fh,offset,ps_sub(jgrid)%w,1,type_subblock_io, &
592 iorequest(iwrite),ierrmpi)
595 if (iwrite>0)
call mpi_waitall(iwrite,iorequest,iostatus,ierrmpi)
596 call mpi_barrier(icomm, ierrmpi)
597 call mpi_file_close(slice_fh,ierrmpi)
600 amode=ior(mpi_mode_append,mpi_mode_wronly)
601 call mpi_file_open(mpi_comm_self,filename,amode,mpi_info_null, &
606 {
call mpi_file_write(slice_fh,subsizes(^de-1),1,mpi_integer,status,ierrmpi)\}
609 call mpi_file_write(slice_fh,nsubleafs,1,mpi_integer,status,ierrmpi)
610 call mpi_file_write(slice_fh,levmax_sub,1,mpi_integer,status,ierrmpi)
611 call mpi_file_write(slice_fh,ndim-1,1,mpi_integer,status,ierrmpi)
612 call mpi_file_write(slice_fh,ndir,1,mpi_integer,status,ierrmpi)
613 call mpi_file_write(slice_fh,nw,1,mpi_integer,status,ierrmpi)
615 call mpi_file_write(slice_fh,it,1,mpi_integer,status,ierrmpi)
616 call mpi_file_write(slice_fh,global_time,1,mpi_double_precision,status,ierrmpi)
618 call mpi_file_close(slice_fh,ierrmpi)
621 end subroutine put_slice_dat
623 subroutine put_slice_zerod
627 integer :: amode, iwrite, status(MPI_STATUS_SIZE)
628 logical,
save :: opened=.false.
629 character(len=name_len) :: wnamei(1:nw+nwauxio),xandwnamei(1:ndim+nw+nwauxio)
630 character(len=1024) :: outfilehead
631 character(len=1024) ::line, data
632 character(len=1024) :: filename, xlabel
633 character(len=79) :: xxlabel
637 write(xlabel,
"(D9.2)")xslice
640 write(xxlabel(1:1),
"(a)")
"+"
642 write(filename,
"(a,i1.1,a,i4.4,a)") trim(base_filename)//
'_d',dir,
'_x'//trim(xxlabel)//
'.csv'
645 if (.not. opened .and. mype==0)
then
646 amode=ior(mpi_mode_create,mpi_mode_wronly)
647 amode=ior(amode,mpi_mode_append)
648 call mpi_file_open(mpi_comm_self,filename,amode, &
649 mpi_info_null,slice_fh,ierrmpi)
653 do iw=1,nw+nwauxio+ndim
654 line = trim(line)//
', '//trim(xandwnamei(iw))
657 call mpi_file_write(slice_fh,line,len_trim(line), &
658 mpi_character,status,ierrmpi)
659 call mpi_file_write(slice_fh,achar(10),1,mpi_character,status,ierrmpi)
660 call mpi_file_close(slice_fh, ierrmpi)
664 call mpi_barrier(icomm,ierrmpi)
668 amode=ior(mpi_mode_create,mpi_mode_wronly)
669 amode=ior(amode,mpi_mode_append)
670 call mpi_file_open(mpi_comm_self,filename,amode, &
671 mpi_info_null,slice_fh,ierrmpi)
674 write(
data,
"(es14.6)")global_time
676 write(
data,
"(es14.6)")ps_sub(1)%x
677 line = trim(line)//
', '//trim(data)
679 write(
data,
"(es14.6)")ps_sub(1)%w(iw)*normconv(iw)
680 line = trim(line)//
', '//trim(data)
683 call mpi_file_write(slice_fh,trim(line),len_trim(line), &
684 mpi_character,status,ierrmpi)
685 call mpi_file_write(slice_fh,achar(10),1,mpi_character,status,ierrmpi)
686 call mpi_file_close(slice_fh, ierrmpi)
689 end subroutine put_slice_zerod
858 integer,
intent(in) :: igrid, dir
859 integer,
intent(inout) :: jgrid
860 logical,
intent(in) :: active
861 double precision,
intent(in) :: xslice
862 double precision,
dimension(0:nw+nwauxio),
intent(out) :: normconv
864 double precision,
dimension(ixMlo^D-1:ixMhi^D,ndim) :: xC_TMP, xC
865 double precision,
dimension(ixMlo^D:ixMhi^D,ndim) :: xCC_TMP, xCC
866 double precision,
dimension(ixMlo^D-1:ixMhi^D,nw+nwauxio) :: wC_TMP
867 double precision,
dimension(ixMlo^D:ixMhi^D,nw+nwauxio) :: wCC_TMP
868 double precision,
dimension(ixG^T) :: x_save
869 integer :: ixslice, nwexpand, ixC^L, ixCC^L
870 logical :: mask(ixG^T)
892 ixslice = minloc(dabs(xslice-ps(igrid)%x(:,dir)),1,mask(:))
897 ixslice = minloc(dabs(xslice-ps(igrid)%x(:,ixmlo2,dir)),1,mask(:,ixmlo2))
899 ixslice = minloc(dabs(xslice-ps(igrid)%x(ixmlo1,:,dir)),1,mask(ixmlo1,:))
901 call mpistop(
"slice direction not clear in fill_subnode")
907 ixslice = minloc(dabs(xslice-ps(igrid)%x(:,ixmlo2,ixmlo3,dir)),1,mask(:,ixmlo2,ixmlo3))
909 ixslice = minloc(dabs(xslice-ps(igrid)%x(ixmlo1,:,ixmlo3,dir)),1,mask(ixmlo1,:,ixmlo3))
911 ixslice = minloc(dabs(xslice-ps(igrid)%x(ixmlo1,ixmlo2,:,dir)),1,mask(ixmlo1,ixmlo2,:))
913 call mpistop(
"slice direction not clear in fill_subnode")
919 xc(:^
d&,dir) = xslice
920 xcc(:^
d&,dir) = xslice
928 ps_sub(jgrid)%w(1:nw+nwexpand) = wcc_tmp(ixslice,1:nw+nwexpand)
929 ps_sub(jgrid)%x(1:
ndim) = xcc_tmp(ixslice,1:
ndim)
930 ps_sub(jgrid)%wC(1:nw+nwexpand) = wc_tmp(ixslice,1:nw+nwexpand)
931 ps_sub(jgrid)%xC(1:
ndim) = xc_tmp(ixslice,1:
ndim)
936 ps_sub(jgrid)%w(ixccmin2:ixccmax2,1:nw+nwexpand) &
937 = wcc_tmp(ixslice,ixccmin2:ixccmax2,1:nw+nwexpand)
938 ps_sub(jgrid)%x(ixccmin2:ixccmax2,1:
ndim) &
939 = xcc_tmp(ixslice,ixccmin2:ixccmax2,1:
ndim)
940 ps_sub(jgrid)%wC(ixcmin2:ixcmax2,1:nw+nwexpand) &
941 = wc_tmp(ixslice,ixcmin2:ixcmax2,1:nw+nwexpand)
942 ps_sub(jgrid)%xC(ixcmin2:ixcmax2,1:
ndim) &
943 = xc_tmp(ixslice,ixcmin2:ixcmax2,1:
ndim)
945 ps_sub(jgrid)%w(ixccmin1:ixccmax1,1:nw+nwexpand) &
946 = wcc_tmp(ixccmin1:ixccmax1,ixslice,1:nw+nwexpand)
947 ps_sub(jgrid)%x(ixccmin1:ixccmax1,1:
ndim) &
948 = xcc_tmp(ixccmin1:ixccmax1,ixslice,1:
ndim)
949 ps_sub(jgrid)%wC(ixcmin1:ixcmax1,1:nw+nwexpand) &
950 = wc_tmp(ixcmin1:ixcmax1,ixslice,1:nw+nwexpand)
951 ps_sub(jgrid)%xC(ixcmin1:ixcmax1,1:
ndim) &
952 = xc_tmp(ixcmin1:ixcmax1,ixslice,1:
ndim)
954 call mpistop(
"slice direction not clear in fill_subnode")
960 ps_sub(jgrid)%w(ixccmin2:ixccmax2,ixccmin3:ixccmax3,1:nw+nwexpand) = &
961 wcc_tmp(ixslice,ixccmin2:ixccmax2,ixccmin3:ixccmax3,1:nw+nwexpand)
962 ps_sub(jgrid)%x(ixccmin2:ixccmax2,ixccmin3:ixccmax3,1:
ndim) = &
963 xcc_tmp(ixslice,ixccmin2:ixccmax2,ixccmin3:ixccmax3,1:
ndim)
964 ps_sub(jgrid)%wC(ixcmin2:ixcmax2,ixcmin3:ixcmax3,1:nw+nwexpand) = &
965 wc_tmp(ixslice,ixcmin2:ixcmax2,ixcmin3:ixcmax3,1:nw+nwexpand)
966 ps_sub(jgrid)%xC(ixcmin2:ixcmax2,ixcmin3:ixcmax3,1:
ndim) = &
967 xc_tmp(ixslice,ixcmin2:ixcmax2,ixcmin3:ixcmax3,1:
ndim)
969 ps_sub(jgrid)%w(ixccmin1:ixccmax1,ixccmin3:ixccmax3,1:nw+nwexpand) = &
970 wcc_tmp(ixccmin1:ixccmax1,ixslice,ixccmin3:ixccmax3,1:nw+nwexpand)
971 ps_sub(jgrid)%x(ixccmin1:ixccmax1,ixccmin3:ixccmax3,1:
ndim) = &
972 xcc_tmp(ixccmin1:ixccmax1,ixslice,ixccmin3:ixccmax3,1:
ndim)
973 ps_sub(jgrid)%wC(ixcmin1:ixcmax1,ixcmin3:ixcmax3,1:nw+nwexpand) = &
974 wc_tmp(ixcmin1:ixcmax1,ixslice,ixcmin3:ixcmax3,1:nw+nwexpand)
975 ps_sub(jgrid)%xC(ixcmin1:ixcmax1,ixcmin3:ixcmax3,1:
ndim) = &
976 xc_tmp(ixcmin1:ixcmax1,ixslice,ixcmin3:ixcmax3,1:
ndim)
978 ps_sub(jgrid)%w(ixccmin1:ixccmax1,ixccmin2:ixccmax2,1:nw+nwexpand) = &
979 wcc_tmp(ixccmin1:ixccmax1,ixccmin2:ixccmax2,ixslice,1:nw+nwexpand)
980 ps_sub(jgrid)%x(ixccmin1:ixccmax1,ixccmin2:ixccmax2,1:
ndim) = &
981 xcc_tmp(ixccmin1:ixccmax1,ixccmin2:ixccmax2,ixslice,1:
ndim)
982 ps_sub(jgrid)%wC(ixcmin1:ixcmax1,ixcmin2:ixcmax2,1:nw+nwexpand) = &
983 wc_tmp(ixcmin1:ixcmax1,ixcmin2:ixcmax2,ixslice,1:nw+nwexpand)
984 ps_sub(jgrid)%xC(ixcmin1:ixcmax1,ixcmin2:ixcmax2,1:
ndim) = &
985 xc_tmp(ixcmin1:ixcmax1,ixcmin2:ixcmax2,ixslice,1:
ndim)
987 call mpistop(
"slice direction not clear in fill_subnode")