diff --git a/src/discretization.f90 b/src/discretization.f90 index 4e60ed54a..3438c6d20 100644 --- a/src/discretization.f90 +++ b/src/discretization.f90 @@ -26,6 +26,9 @@ module discretization discretization_NodeCoords0, & discretization_IPcoords, & discretization_NodeCoords + + integer :: & + discretization_sharedNodesBegin public :: & discretization_init, & @@ -38,7 +41,9 @@ contains !-------------------------------------------------------------------------------------------------- !> @brief stores the relevant information in globally accesible variables !-------------------------------------------------------------------------------------------------- -subroutine discretization_init(homogenizationAt,microstructureAt,IPcoords0,NodeCoords0) +subroutine discretization_init(homogenizationAt,microstructureAt,& + IPcoords0,NodeCoords0,& + sharedNodesBegin) integer, dimension(:), intent(in) :: & homogenizationAt, & @@ -46,6 +51,8 @@ subroutine discretization_init(homogenizationAt,microstructureAt,IPcoords0,NodeC real(pReal), dimension(:,:), intent(in) :: & IPcoords0, & NodeCoords0 + integer, optional, intent(in) :: & + sharedNodesBegin write(6,'(/,a)') ' <<<+- discretization init -+>>>' @@ -61,6 +68,12 @@ subroutine discretization_init(homogenizationAt,microstructureAt,IPcoords0,NodeC discretization_NodeCoords0 = NodeCoords0 discretization_NodeCoords = NodeCoords0 + if(present(sharedNodesBegin)) then + discretization_sharedNodesBegin = sharedNodesBegin + else + discretization_sharedNodesBegin = size(discretization_NodeCoords0,2) + endif + end subroutine discretization_init @@ -73,10 +86,12 @@ subroutine discretization_results call HDF5_closeGroup(results_addGroup(trim('current/geometry'))) - u = discretization_NodeCoords - discretization_NodeCoords0 + u = discretization_NodeCoords (1:3,:discretization_sharedNodesBegin) & + - discretization_NodeCoords0(1:3,:discretization_sharedNodesBegin) call results_writeDataset('current/geometry',u,'u_n','nodal displacements','m') - u = discretization_IPcoords - discretization_IPcoords0 + u = discretization_IPcoords & + - discretization_IPcoords0 call results_writeDataset('current/geometry',u,'u_c','cell center displacements','m') #endif end subroutine discretization_results diff --git a/src/mesh_grid.f90 b/src/mesh_grid.f90 index a69c9b807..ca25e34d4 100644 --- a/src/mesh_grid.f90 +++ b/src/mesh_grid.f90 @@ -91,7 +91,10 @@ subroutine mesh_init(ip,el) call discretization_init(homogenizationAt,microstructureAt, & IPcoordinates0(myGrid,mySize,grid3Offset), & - Nodes0(myGrid,mySize,grid3Offset)) + Nodes0(myGrid,mySize,grid3Offset),& + merge((grid(1)+1) * (grid(2)+1) * (grid3+1),& ! write bottom layer + (grid(1)+1) * (grid(2)+1) * grid3,& ! do not write bottom layer (is top of rank-1) + worldrank<1)) FEsolving_execElem = [1,product(myGrid)] ! parallel loop bounds set to comprise all elements allocate(FEsolving_execIP(2,product(myGrid)),source=1) ! parallel loop bounds set to comprise the only IP