xref: /phasta/phSolver/incompressible/genlmass.f (revision 595995161822a203c8467e0e4a253d7bd7d6df32)
1*59599516SKenneth E. Jansen        subroutine genlmass (x, shp,shgl)
2*59599516SKenneth E. Jansenc
3*59599516SKenneth E. Jansen        use pointer_data
4*59599516SKenneth E. Jansenc
5*59599516SKenneth E. Jansen        include "common.h"
6*59599516SKenneth E. Jansen        include "mpif.h"
7*59599516SKenneth E. Jansenc
8*59599516SKenneth E. Jansen        real*8 x(numnp,nsd)
9*59599516SKenneth E. Jansenc
10*59599516SKenneth E. Jansen        real*8 shp(MAXTOP,maxsh,MAXQPT),
11*59599516SKenneth E. Jansen     &            shgl(MAXTOP,nsd,maxsh,MAXQPT)
12*59599516SKenneth E. Jansenc
13*59599516SKenneth E. Jansen        real*8, allocatable :: tmpshp(:,:), tmpshgl(:,:,:)
14*59599516SKenneth E. Jansen
15*59599516SKenneth E. Jansenc
16*59599516SKenneth E. Jansenc gmass came in via pointer_data and will
17*59599516SKenneth E. Jansenc be available wherever it is included  (allocate it now).
18*59599516SKenneth E. Jansenc
19*59599516SKenneth E. Jansen
20*59599516SKenneth E. Jansen        allocate (gmass(nshg))
21*59599516SKenneth E. Jansen        gmass=zero
22*59599516SKenneth E. Jansenc
23*59599516SKenneth E. Jansenc.... loop over the element-blocks
24*59599516SKenneth E. Jansenc
25*59599516SKenneth E. Jansen        do iblk = 1, nelblk
26*59599516SKenneth E. Jansen          iel    = lcblk(1,iblk)
27*59599516SKenneth E. Jansen          lelCat = lcblk(2,iblk)
28*59599516SKenneth E. Jansen          lcsyst = lcblk(3,iblk)
29*59599516SKenneth E. Jansen          iorder = lcblk(4,iblk)
30*59599516SKenneth E. Jansen          nenl   = lcblk(5,iblk) ! no. of vertices per element
31*59599516SKenneth E. Jansen          nshl   = lcblk(10,iblk)
32*59599516SKenneth E. Jansen          mattyp = lcblk(7,iblk)
33*59599516SKenneth E. Jansen          npro   = lcblk(1,iblk+1) - iel
34*59599516SKenneth E. Jansen          inum   = iel + npro - 1
35*59599516SKenneth E. Jansen          ngauss = nint(lcsyst)
36*59599516SKenneth E. Jansenc
37*59599516SKenneth E. Jansenc
38*59599516SKenneth E. Jansenc.... compute and assemble the residual and tangent matrix
39*59599516SKenneth E. Jansenc
40*59599516SKenneth E. Jansen          allocate (tmpshp(nshl,MAXQPT))
41*59599516SKenneth E. Jansen          allocate (tmpshgl(nsd,nshl,MAXQPT))
42*59599516SKenneth E. Jansen          tmpshp(1:nshl,:) = shp(lcsyst,1:nshl,:)
43*59599516SKenneth E. Jansen          tmpshgl(:,1:nshl,:) = shgl(lcsyst,:,1:nshl,:)
44*59599516SKenneth E. Jansen
45*59599516SKenneth E. Jansen
46*59599516SKenneth E. Jansen          call AsImass (x,       tmpshp,
47*59599516SKenneth E. Jansen     &                  tmpshgl, mien(iblk)%p,
48*59599516SKenneth E. Jansen     &                  gmass)
49*59599516SKenneth E. Jansen
50*59599516SKenneth E. Jansen          deallocate ( tmpshp )
51*59599516SKenneth E. Jansen          deallocate ( tmpshgl )
52*59599516SKenneth E. Jansenc
53*59599516SKenneth E. Jansenc.... end of interior element loop
54*59599516SKenneth E. Jansenc
55*59599516SKenneth E. Jansen       enddo
56*59599516SKenneth E. Jansen
57*59599516SKenneth E. Jansen      return
58*59599516SKenneth E. Jansen      end
59*59599516SKenneth E. Jansen
60*59599516SKenneth E. Jansen        subroutine AsImass (x,      shp,
61*59599516SKenneth E. Jansen     &                     shgl,    ien,
62*59599516SKenneth E. Jansen     &                     gmass)
63*59599516SKenneth E. Jansenc
64*59599516SKenneth E. Jansenc----------------------------------------------------------------------
65*59599516SKenneth E. Jansenc
66*59599516SKenneth E. Jansenc This routine computes and assembles the mass corresponding to the
67*59599516SKenneth E. Jansenc each node.
68*59599516SKenneth E. Jansenc
69*59599516SKenneth E. Jansenc Ken Jansen, Winter 2000.  (Fortran 90)
70*59599516SKenneth E. Jansenc----------------------------------------------------------------------
71*59599516SKenneth E. Jansenc
72*59599516SKenneth E. Jansen      include "common.h"
73*59599516SKenneth E. Jansenc
74*59599516SKenneth E. Jansen        real*8 x(numnp,nsd),
75*59599516SKenneth E. Jansen     &         shp(nshl,maxsh),       shgl(nsd,nshl,ngauss),
76*59599516SKenneth E. Jansen     &         gmass(nshg)
77*59599516SKenneth E. Jansen
78*59599516SKenneth E. Jansen        integer ien(npro,nshl)
79*59599516SKenneth E. Jansen
80*59599516SKenneth E. Jansenc
81*59599516SKenneth E. Jansen        real*8    xl(npro,nenl,nsd),    WdetJ(npro),
82*59599516SKenneth E. Jansen     &            sgn(npro,nshl),       shape(npro,nshl),
83*59599516SKenneth E. Jansen     &            locmass(npro,nshl),   shg(npro,nshl,nsd),
84*59599516SKenneth E. Jansen     &            fmstot(npro),         temp(npro),
85*59599516SKenneth E. Jansen     &            dxidx(npro,nsd,nsd),  shdrv(npro,nsd,nshl)
86*59599516SKenneth E. Jansen
87*59599516SKenneth E. Jansen        integer aa
88*59599516SKenneth E. Jansenc
89*59599516SKenneth E. Jansenc
90*59599516SKenneth E. Jansenc
91*59599516SKenneth E. Jansenc.... gather the variables
92*59599516SKenneth E. Jansenc
93*59599516SKenneth E. Jansenc
94*59599516SKenneth E. Jansenc.... get the matrix of mode signs for the hierarchic basis functions.
95*59599516SKenneth E. Jansenc
96*59599516SKenneth E. Jansen        if (ipord .gt. 1) then
97*59599516SKenneth E. Jansen           call getsgn(ien,sgn)
98*59599516SKenneth E. Jansen        endif
99*59599516SKenneth E. Jansen
100*59599516SKenneth E. Jansen        call localx(x,      xl,     ien,    nsd,    'gather  ')
101*59599516SKenneth E. Jansenc
102*59599516SKenneth E. Jansenc.... zero the matrices if they are being recalculated
103*59599516SKenneth E. Jansenc
104*59599516SKenneth E. Jansen
105*59599516SKenneth E. Jansen        locmass=zero
106*59599516SKenneth E. Jansen        fmstot=zero
107*59599516SKenneth E. Jansen
108*59599516SKenneth E. Jansen        do intp = 1, ngauss
109*59599516SKenneth E. Jansen
110*59599516SKenneth E. Jansen           if (Qwt(lcsyst,intp) .eq. zero) cycle ! precaution
111*59599516SKenneth E. Jansenc
112*59599516SKenneth E. Jansenc.... get the hierarchic shape functions at this int point
113*59599516SKenneth E. Jansenc
114*59599516SKenneth E. Jansen           call getshp(shp,          shgl,      sgn,
115*59599516SKenneth E. Jansen     &                 shape,        shdrv)
116*59599516SKenneth E. Jansen
117*59599516SKenneth E. Jansenc
118*59599516SKenneth E. Jansenc.... --------------------->  Element Metrics  <-----------------------
119*59599516SKenneth E. Jansenc
120*59599516SKenneth E. Jansen           call e3metric( xl,         shdrv,       dxidx,
121*59599516SKenneth E. Jansen     &                    shg,        WdetJ)
122*59599516SKenneth E. Jansen
123*59599516SKenneth E. Jansenc
124*59599516SKenneth E. Jansenc  get this quad points contribution to the integral of the square of  the
125*59599516SKenneth E. Jansenc  shape function
126*59599516SKenneth E. Jansenc
127*59599516SKenneth E. Jansen           do aa = 1,nshl
128*59599516SKenneth E. Jansen              locmass(:,aa)= locmass(:,aa)
129*59599516SKenneth E. Jansen     &             + shape(:,aa)*shape(:,aa)*WdetJ
130*59599516SKenneth E. Jansen           enddo
131*59599516SKenneth E. Jansenc
132*59599516SKenneth E. Jansenc also accumulate this quad points contribution to the integral of the element
133*59599516SKenneth E. Jansenc volume (integral Na^2 d Omega)
134*59599516SKenneth E. Jansenc
135*59599516SKenneth E. Jansen           fmstot= fmstot + WdetJ ! intregral  d Omega
136*59599516SKenneth E. Jansenc
137*59599516SKenneth E. Jansenc.... end of integration loop
138*59599516SKenneth E. Jansenc
139*59599516SKenneth E. Jansen        enddo
140*59599516SKenneth E. Jansenc
141*59599516SKenneth E. Jansenc.... lumped mass if needed   Note that the locmass factors accumulated
142*59599516SKenneth E. Jansenc     over integration points and weighted with WdetJ already.
143*59599516SKenneth E. Jansenc
144*59599516SKenneth E. Jansen
145*59599516SKenneth E. Jansenc.... scale the LHS matrix contribution with special lumping weighting
146*59599516SKenneth E. Jansenc
147*59599516SKenneth E. Jansenc  The first term we collect is the trace of integral Na^2 d Omega
148*59599516SKenneth E. Jansenc
149*59599516SKenneth E. Jansen        temp = zero
150*59599516SKenneth E. Jansen        do aa = 1, nshl
151*59599516SKenneth E. Jansen           temp = temp + locmass(:,aa) !reusing temp to save memory
152*59599516SKenneth E. Jansen        enddo
153*59599516SKenneth E. Jansen
154*59599516SKenneth E. Jansenc
155*59599516SKenneth E. Jansenc scale the diagonal so that the trace will still yield Omega^e (the volume
156*59599516SKenneth E. Jansenc of the element)
157*59599516SKenneth E. Jansenc
158*59599516SKenneth E. Jansen        do aa = 1, nshl
159*59599516SKenneth E. Jansen           locmass(:,aa) = locmass(:,aa) * fmstot / temp
160*59599516SKenneth E. Jansen        enddo
161*59599516SKenneth E. Jansenc
162*59599516SKenneth E. Jansenc.... assemble the residual
163*59599516SKenneth E. Jansenc
164*59599516SKenneth E. Jansen        call local (gmass,    locmass,     ien,    1,  'scatter ')
165*59599516SKenneth E. Jansen
166*59599516SKenneth E. Jansenc
167*59599516SKenneth E. Jansenc.... end
168*59599516SKenneth E. Jansenc
169*59599516SKenneth E. Jansen        return
170*59599516SKenneth E. Jansen        end
171*59599516SKenneth E. Jansen
172*59599516SKenneth E. Jansen      subroutine lmassadd ( ac,       res,
173*59599516SKenneth E. Jansen     &                      rowp,     colm,
174*59599516SKenneth E. Jansen     &                      lhsK,     gmass)
175*59599516SKenneth E. Jansenc
176*59599516SKenneth E. Jansen      include "common.h"
177*59599516SKenneth E. Jansenc
178*59599516SKenneth E. Jansen      real*8 ac(nshg,ndof), res(nshg,4), tmp,tmp1
179*59599516SKenneth E. Jansen      real*8 lhsK(9,nnz_tot), gmass(nshg), rho(nshg)
180*59599516SKenneth E. Jansen      integer rowp(nnz*nshg),  colm(nshg+1)
181*59599516SKenneth E. Jansen      integer	n,	k
182*59599516SKenneth E. Jansenc
183*59599516SKenneth E. Jansen      integer sparseloc
184*59599516SKenneth E. Jansenc
185*59599516SKenneth E. Jansenc
186*59599516SKenneth E. Jansen      rho=datmat(1,1,1)  ! needs to be generalized for VOF or level set
187*59599516SKenneth E. Jansen      tmp1=flmpl*almi
188*59599516SKenneth E. Jansen      if((flmpl.ne.0).and.(lhs.eq.1)) then
189*59599516SKenneth E. Jansenc
190*59599516SKenneth E. Jansenc.... Add lmass to diag of lhsK
191*59599516SKenneth E. Jansenc
192*59599516SKenneth E. Jansen         do n = 1, nshg
193*59599516SKenneth E. Jansen	    k = sparseloc( rowp(colm(n)), colm(n+1)-colm(n), n )
194*59599516SKenneth E. Jansen     &       + colm(n)-1
195*59599516SKenneth E. Jansen            tmp=gmass(n)*tmp1*rho(n)
196*59599516SKenneth E. Jansen	    lhsK(1,k) = lhsK(1,k) + tmp
197*59599516SKenneth E. Jansen	    lhsK(5,k) = lhsK(5,k) + tmp
198*59599516SKenneth E. Jansen	    lhsK(9,k) = lhsK(9,k) + tmp
199*59599516SKenneth E. Jansen         enddo
200*59599516SKenneth E. Jansen      endif
201*59599516SKenneth E. Jansen
202*59599516SKenneth E. Jansen      tmp1=flmpr
203*59599516SKenneth E. Jansen
204*59599516SKenneth E. Jansen      if(flmpr.ne.0) then
205*59599516SKenneth E. Jansen         rho=rho*gmass*tmp1  ! reuse rho
206*59599516SKenneth E. Jansen         res(:,1)=res(:,1)-ac(:,1)*rho(:)
207*59599516SKenneth E. Jansen         res(:,2)=res(:,2)-ac(:,2)*rho(:)
208*59599516SKenneth E. Jansen         res(:,3)=res(:,3)-ac(:,3)*rho(:)
209*59599516SKenneth E. Jansen      endif
210*59599516SKenneth E. Jansen
211*59599516SKenneth E. Jansen      return
212*59599516SKenneth E. Jansen      end
213*59599516SKenneth E. Jansen
214*59599516SKenneth E. Jansen      subroutine lmassaddSclr ( ac,       res,
215*59599516SKenneth E. Jansen     &                          rowp,     colm,
216*59599516SKenneth E. Jansen     &                          lhsS,     gmass)
217*59599516SKenneth E. Jansenc
218*59599516SKenneth E. Jansen      include "common.h"
219*59599516SKenneth E. Jansenc
220*59599516SKenneth E. Jansen      real*8 ac(nshg),       res(nshg), tmp, tmp1
221*59599516SKenneth E. Jansen      real*8 lhsS(nnz_tot), gmass(nshg), rho(nshg)
222*59599516SKenneth E. Jansen      integer rowp(nnz*nshg),  colm(nshg+1)
223*59599516SKenneth E. Jansen      integer	n,	k
224*59599516SKenneth E. Jansenc
225*59599516SKenneth E. Jansen      integer sparseloc
226*59599516SKenneth E. Jansenc
227*59599516SKenneth E. Jansenc
228*59599516SKenneth E. Jansen      rho=datmat(1,1,1)  ! needs to be generalized for VOF or level set
229*59599516SKenneth E. Jansen      tmp1=flmpl*almi
230*59599516SKenneth E. Jansen      if((flmpl.ne.0).and.(lhs.eq.1)) then
231*59599516SKenneth E. Jansenc
232*59599516SKenneth E. Jansenc.... Add lmass to diag of lhsK
233*59599516SKenneth E. Jansenc
234*59599516SKenneth E. Jansen         do n = 1, nshg
235*59599516SKenneth E. Jansen	    k = sparseloc( rowp(colm(n)), colm(n+1)-colm(n), n )
236*59599516SKenneth E. Jansen     &       + colm(n)-1
237*59599516SKenneth E. Jansen            tmp=gmass(n)*tmp1*rho(n)
238*59599516SKenneth E. Jansen	    lhsS(k) = lhsS(k) + tmp
239*59599516SKenneth E. Jansen         enddo
240*59599516SKenneth E. Jansen      endif
241*59599516SKenneth E. Jansen
242*59599516SKenneth E. Jansen      tmp1=flmpr
243*59599516SKenneth E. Jansen      if(flmpr.ne.0) then
244*59599516SKenneth E. Jansen         rho=rho*gmass*tmp1  ! reuse rho
245*59599516SKenneth E. Jansen         res(:)=res(:)-ac(:)*rho(:)
246*59599516SKenneth E. Jansen      endif
247*59599516SKenneth E. Jansen
248*59599516SKenneth E. Jansen      return
249*59599516SKenneth E. Jansen      end
250