File: resultsforc.f

package info (click to toggle)
calculix-ccx 2.11-1
  • links: PTS, VCS
  • area: main
  • in suites: buster, stretch
  • size: 10,188 kB
  • sloc: fortran: 115,312; ansic: 34,480; sh: 374; makefile: 35; perl: 15
file content (118 lines) | stat: -rw-r--r-- 4,014 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
!
!     CalculiX - A 3-dimensional finite element program
!              Copyright (C) 1998-2015 Guido Dhondt
!
!     This program is free software; you can redistribute it and/or
!     modify it under the terms of the GNU General Public License as
!     published by the Free Software Foundation(version 2);
!     
!
!     This program is distributed in the hope that it will be useful,
!     but WITHOUT ANY WARRANTY; without even the implied warranty of 
!     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 
!     GNU General Public License for more details.
!
!     You should have received a copy of the GNU General Public License
!     along with this program; if not, write to the Free Software
!     Foundation, Inc., 675 Mass Ave, Cambridge, MA 02139, USA.
!
      subroutine resultsforc(nk,f,fn,nactdof,ipompc,nodempc,
     &  coefmpc,labmpc,nmpc,mi,fmpc,calcul_fn,calcul_f)
!
!     calculating the equation system internal force vector
!     (one entry for each active degree of freedom)
!
      implicit none
!
      character*20 labmpc(*)
!
      integer mi(*),nactdof(0:mi(2),*),ipompc(*),nodempc(3,*),nk,i,j,
     &  nmpc,ist,ndir,node,index,calcul_fn,calcul_f
!
      real*8 f(*),fn(0:mi(2),*),coefmpc(*),fmpc(*),forcempc
!
!     subtracting the mpc force (for each linear mpc there is one
!     force; the actual force in a node belonging to the mpc is
!     obtained by multiplying this force with the nodal coefficient.
!     The force has to be subtracted from f, since it does not
!     appear on the rhs of the equation system
!
      if(calcul_fn.eq.1)then
        do i=1,nmpc
            ist=ipompc(i)
            node=nodempc(1,ist)
            ndir=nodempc(2,ist)
            if(ndir.gt.3) cycle
            forcempc=fn(ndir,node)/coefmpc(ist)
            fmpc(i)=forcempc
            fn(ndir,node)=0.d0
            index=nodempc(3,ist)
            if(index.eq.0) cycle
            do
               node=nodempc(1,index)
               ndir=nodempc(2,index)
               fn(ndir,node)=fn(ndir,node)-coefmpc(index)*forcempc
               index=nodempc(3,index)
               if(index.eq.0) exit
            enddo
         enddo
      endif
!
!     calculating the system force vector
!
      if(calcul_f.eq.1) then
         do i=1,nk
            do j=0,mi(2)
               if(nactdof(j,i).gt.0) then
                  f(nactdof(j,i))=fn(j,i)
               endif
            enddo
         enddo
      endif
!
!     adding the mpc force again to fn
!
      if(calcul_fn.eq.1)then
         do i=1,nmpc
            ist=ipompc(i)
            node=nodempc(1,ist)
            ndir=nodempc(2,ist)
            if(ndir.gt.3) cycle
            forcempc=fmpc(i)
            fn(ndir,node)=forcempc*coefmpc(ist)
            index=nodempc(3,ist)
!
!           nodes not belonging to the structure have to be
!           taken out
!
            if(labmpc(i)(1:7).eq.'MEANROT') then
               if(nodempc(3,nodempc(3,index)).eq.0) cycle
            elseif(labmpc(i)(1:10).eq.'PRETENSION') then
               if(nodempc(3,index).eq.0) cycle
            elseif(labmpc(i)(1:5).eq.'RIGID') then
               if(nodempc(3,nodempc(3,nodempc(3,nodempc(3,nodempc(3,inde
     &x))))).eq.0) cycle
            else
               if(index.eq.0) cycle
            endif
            do
               node=nodempc(1,index)
               ndir=nodempc(2,index)
               fn(ndir,node)=fn(ndir,node)+coefmpc(index)*forcempc
               index=nodempc(3,index)
               if(labmpc(i)(1:7).eq.'MEANROT') then
                  if(nodempc(3,nodempc(3,index)).eq.0) exit
               elseif(labmpc(i)(1:10).eq.'PRETENSION') then
                  if(nodempc(3,index).eq.0) exit
               elseif(labmpc(i)(1:5).eq.'RIGID') then
                  if(nodempc(3,nodempc(3,nodempc(3,nodempc(3,nodempc(3,i
     &ndex))))).eq.0) exit
               else
                  if(index.eq.0) exit
               endif
            enddo
         enddo
      endif
!
      return
      end