little polishing
- same names as for grid - 'local' naming misleading here (probably for x_local etc also not a good choice) - following names in PETSc notation (https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMPLEX/DMPlexGetPointLocal.html) - using allocation with bounds to simplify access - NOTE: The nodal coordinates are not what damask.Result expects except for simple, 1 IP elements!
This commit is contained in:
parent
1516ded71b
commit
a0ffbc5f17
|
@ -372,7 +372,7 @@ program DAMASK_mesh
|
||||||
|
|
||||||
if (mod(inc,loadCases(currentLoadCase)%outputFrequency) == 0) then ! at output frequency
|
if (mod(inc,loadCases(currentLoadCase)%outputFrequency) == 0) then ! at output frequency
|
||||||
print'(/,a)', ' ... writing results to file ......................................'
|
print'(/,a)', ' ... writing results to file ......................................'
|
||||||
call FEM_mechanical_updatenodeandipcoords
|
call FEM_mechanical_updateCoords
|
||||||
call CPFEM_results(totalIncsCounter,time)
|
call CPFEM_results(totalIncsCounter,time)
|
||||||
endif
|
endif
|
||||||
|
|
||||||
|
|
|
@ -69,7 +69,7 @@ module mesh_mechanical_FEM
|
||||||
FEM_mechanical_init, &
|
FEM_mechanical_init, &
|
||||||
FEM_mechanical_solution, &
|
FEM_mechanical_solution, &
|
||||||
FEM_mechanical_forward, &
|
FEM_mechanical_forward, &
|
||||||
FEM_mechanical_updatenodeandipcoords
|
FEM_mechanical_updateCoords
|
||||||
|
|
||||||
contains
|
contains
|
||||||
|
|
||||||
|
@ -666,43 +666,40 @@ subroutine FEM_mechanical_converged(snes_local,PETScIter,xnorm,snorm,fnorm,reaso
|
||||||
|
|
||||||
end subroutine FEM_mechanical_converged
|
end subroutine FEM_mechanical_converged
|
||||||
|
|
||||||
|
|
||||||
!--------------------------------------------------------------------------------------------------
|
!--------------------------------------------------------------------------------------------------
|
||||||
!> @brief write diplacements
|
!> @brief Calculate current coordinates (FEM nodal coordinates only at the moment)
|
||||||
!--------------------------------------------------------------------------------------------------
|
!--------------------------------------------------------------------------------------------------
|
||||||
subroutine FEM_mechanical_updatenodeandipcoords
|
subroutine FEM_mechanical_updateCoords()
|
||||||
|
|
||||||
real(pReal), pointer, dimension(:) :: &
|
real(pReal), pointer, dimension(:) :: &
|
||||||
node_coords_local !nodal coordinates (shape - (dimplex,Nnodes))
|
nodeCoords_linear !< nodal coordinates (dimPlex*Nnodes)
|
||||||
|
|
||||||
real(pReal), pointer, dimension(:,:) :: &
|
real(pReal), pointer, dimension(:,:) :: &
|
||||||
node_coords !nodal coordinates (shape - (3,Nnodes), i.e. with dummy values for 2D elements)
|
nodeCoords !< nodal coordinates (3,Nnodes)
|
||||||
|
|
||||||
DM :: dm_local
|
DM :: dm_local
|
||||||
Vec :: x_local
|
Vec :: x_local
|
||||||
PetscErrorCode :: ierr
|
PetscErrorCode :: ierr
|
||||||
PetscInt :: dimPlex, vStart, vEnd, v, x, y
|
PetscInt :: dimPlex, pStart, pEnd, p, s, e
|
||||||
PetscInt :: v_count = 1
|
|
||||||
PetscSection :: section
|
PetscSection :: section
|
||||||
|
|
||||||
call SNESGetDM(mechanical_snes,dm_local,ierr); CHKERRQ(ierr)
|
call SNESGetDM(mechanical_snes,dm_local,ierr); CHKERRQ(ierr)
|
||||||
call DMGetLocalSection(dm_local,section,ierr); CHKERRQ(ierr)
|
call DMGetLocalSection(dm_local,section,ierr); CHKERRQ(ierr)
|
||||||
call DMGetLocalVector(dm_local,x_local,ierr); CHKERRQ(ierr)
|
call DMGetLocalVector(dm_local,x_local,ierr); CHKERRQ(ierr)
|
||||||
call DMGetDimension(dm_local,dimPlex,ierr); CHKERRQ(ierr)
|
call DMGetDimension(dm_local,dimPlex,ierr); CHKERRQ(ierr)
|
||||||
call DMPlexGetDepthStratum(dm_local, 0, vStart, vEnd,ierr); CHKERRQ(ierr)
|
call DMPlexGetDepthStratum(dm_local,0,pStart,pEnd,ierr); CHKERRQ(ierr)
|
||||||
allocate(node_coords(3,vEnd-vStart),source=0.0_pReal)
|
allocate(nodeCoords(3,pStart:pEnd-1),source=0.0_pReal)
|
||||||
call VecGetArrayF90(x_local, node_coords_local,ierr); CHKERRQ(ierr)
|
call VecGetArrayF90(x_local,nodeCoords_linear,ierr); CHKERRQ(ierr)
|
||||||
|
|
||||||
v_count=1
|
do p=pStart, pEnd-1
|
||||||
do v=vStart, vEnd-1 !Loop over vertices
|
call DMPlexGetPointLocal(dm_local, p, s, e, ierr); CHKERRQ(ierr)
|
||||||
call DMPlexGetPointLocal(dm_local, v, x, y, ierr); CHKERRQ(ierr)
|
nodeCoords(1:dimPlex,p)=nodeCoords_linear(s+1:e)
|
||||||
node_coords(1:dimPlex,v_count)=node_coords_local(x+1:y)
|
|
||||||
v_count = v_count+1
|
|
||||||
end do
|
end do
|
||||||
|
|
||||||
call discretization_setNodeCoords(node_coords)
|
call discretization_setNodeCoords(nodeCoords)
|
||||||
call VecRestoreArrayF90(x_local,node_coords_local,ierr); CHKERRQ(ierr)
|
call VecRestoreArrayF90(x_local,nodeCoords_linear,ierr); CHKERRQ(ierr)
|
||||||
call DMRestoreLocalVector(dm_local,x_local,ierr); CHKERRQ(ierr)
|
call DMRestoreLocalVector(dm_local,x_local,ierr); CHKERRQ(ierr)
|
||||||
|
|
||||||
end subroutine FEM_mechanical_updatenodeandipcoords
|
end subroutine FEM_mechanical_updateCoords
|
||||||
|
|
||||||
end module mesh_mechanical_FEM
|
end module mesh_mechanical_FEM
|
||||||
|
|
Loading…
Reference in New Issue