FrontISTR 5.2.0
Large-scale structural analysis program with finit element method
Loading...
Searching...
No Matches
fstr_ctrl_modifier.f90
Go to the documentation of this file.
1!-------------------------------------------------------------------------------
2! Copyright (c) 2019 FrontISTR Commons
3! This software is released under the MIT License, see LICENSE.txt
4!-------------------------------------------------------------------------------
7
9 use hecmw_util
10 implicit none
11
12contains
13
15 subroutine fstr_append_mpc( np, nodes, dofs, values, mpcs )
16 integer, intent(in) :: np
17 integer, intent(in) :: nodes(np)
18 integer, intent(in) :: dofs(np)
19 real(kind=kreal), intent(in) :: values(np+1)
20 type( hecmwst_mpc ), intent(inout) :: mpcs
21
22 integer :: i, n_mpc, old_size, new_size
23 n_mpc = mpcs%n_mpc
24 new_size= n_mpc+1
25 mpcs%n_mpc = new_size
26 call fstr_expand_index_array( mpcs%mpc_index, n_mpc+1, new_size+1 )
27 call fstr_expand_real_array( mpcs%mpc_const, n_mpc, n_mpc+1 )
28 old_size = mpcs%mpc_index( n_mpc )
29 new_size = old_size+np
30 call fstr_expand_integer_array( mpcs%mpc_item, old_size, new_size )
31 call fstr_expand_integer_array( mpcs%mpc_dof, old_size, new_size )
32 call fstr_expand_real_array( mpcs%mpc_val, old_size, new_size )
33
34 mpcs%mpc_index(mpcs%n_mpc) = mpcs%mpc_index(mpcs%n_mpc-1)+np
35 mpcs%mpc_const(mpcs%n_mpc) = values(np+1)
36 do i=1,np
37 mpcs%mpc_item(old_size+i) = nodes(i)
38 mpcs%mpc_dof(old_size+i) = dofs(i)
39 mpcs%mpc_val(old_size+i) = values(i)
40 enddo
41 end subroutine
42
44 subroutine fstr_delete_mpc( np, mpcs )
45 integer, intent(in) :: np
46 type( hecmwst_mpc ), intent(inout) :: mpcs
47
48 integer :: n_mpc, old_size, nitem
49 n_mpc = mpcs%n_mpc
50 old_size = mpcs%mpc_index( n_mpc )
51 nitem = old_size- mpcs%mpc_index( n_mpc-np )
52 call fstr_delete_real_array( mpcs%mpc_val, old_size, nitem )
53 call fstr_delete_integer_array( mpcs%mpc_dof, old_size, nitem )
54 call fstr_delete_integer_array( mpcs%mpc_item, old_size, nitem )
55
56 call fstr_delete_real_array( mpcs%mpc_const, n_mpc, np )
57 call fstr_delete_index_array( mpcs%mpc_index, n_mpc, np )
58 mpcs%n_mpc = n_mpc-np
59 end subroutine
60
61end module
62
This module provides functions to modify MPC conditions.
subroutine fstr_delete_mpc(np, mpcs)
Delete last n equation conditions from current mpc condition.
subroutine fstr_append_mpc(np, nodes, dofs, values, mpcs)
Append new equation condition at end of existing mpc conditions.
This module contains auxiliary functions in calculation setup.
subroutine fstr_delete_real_array(array, old_size, nitem)
subroutine fstr_expand_integer_array(array, old_size, new_size)
subroutine fstr_expand_index_array(array, old_size, new_size)
subroutine fstr_expand_real_array(array, old_size, new_size)
subroutine fstr_delete_integer_array(array, old_size, nitem)
subroutine fstr_delete_index_array(array, old_size, nindex)
I/O and Utility.
Definition: hecmw_util_f.F90:7
integer(kind=4), parameter kreal