1*59599516SKenneth E. Jansen subroutine solvecon(y, x, iBC, BC, 2*59599516SKenneth E. Jansen & iper, ilwork, shp, shgl) 3*59599516SKenneth E. Jansen 4*59599516SKenneth E. Jansenc--------------------------------------------------------------------- 5*59599516SKenneth E. Jansenc This subroutine is to calculate the constarint for the redistancing 6*59599516SKenneth E. Jansenc scalar of the level set method. This is to prevent interface from 7*59599516SKenneth E. Jansenc moving by applying the condition that the volume must stay constant 8*59599516SKenneth E. Jansenc in each element when the redisatnce step is applied. 9*59599516SKenneth E. Jansenc-------------------------------------------------------------------- 10*59599516SKenneth E. Jansenc 11*59599516SKenneth E. Jansenc 12*59599516SKenneth E. Jansen use pointer_data 13*59599516SKenneth E. Jansenc 14*59599516SKenneth E. Jansen include "common.h" 15*59599516SKenneth E. Jansen include "mpif.h" 16*59599516SKenneth E. Jansen include "auxmpi.h" 17*59599516SKenneth E. Jansenc 18*59599516SKenneth E. Jansen dimension y(nshg,ndof), 19*59599516SKenneth E. Jansen & x(numnp,nsd), iBC(nshg), 20*59599516SKenneth E. Jansen & BC(nshg,ndofBC), ilwork(nlwork), 21*59599516SKenneth E. Jansen & iper(nshg) 22*59599516SKenneth E. Jansenc 23*59599516SKenneth E. Jansenc 24*59599516SKenneth E. Jansen dimension shp(MAXTOP,maxsh,MAXQPT), 25*59599516SKenneth E. Jansen & shgl(MAXTOP,nsd,maxsh,MAXQPT), 26*59599516SKenneth E. Jansen & v_lambda(nshg), hprime(nshg), 27*59599516SKenneth E. Jansen & v_lambda1(nshg), v_lambda2(nshg), 28*59599516SKenneth E. Jansen & rmass(nshg) 29*59599516SKenneth E. Jansenc 30*59599516SKenneth E. Jansen real*8, allocatable :: tmpshp(:,:), tmpshgl(:,:,:) 31*59599516SKenneth E. Jansen real*8, allocatable :: tmpshpb(:,:), tmpshglb(:,:,:) 32*59599516SKenneth E. Jansen 33*59599516SKenneth E. Jansenc 34*59599516SKenneth E. Jansenc ... intialize 35*59599516SKenneth E. Jansenc 36*59599516SKenneth E. Jansen rmass = zero 37*59599516SKenneth E. Jansen v_lambda = zero 38*59599516SKenneth E. Jansen v_lambda1 = zero 39*59599516SKenneth E. Jansen v_lambda2 = zero 40*59599516SKenneth E. Jansen hprime = zero 41*59599516SKenneth E. Jansenc 42*59599516SKenneth E. Jansenc ... loop over element blocks 43*59599516SKenneth E. Jansenc 44*59599516SKenneth E. Jansen do iblk = 1, nelblk 45*59599516SKenneth E. Jansenc 46*59599516SKenneth E. Jansenc.... set up the parameters 47*59599516SKenneth E. Jansenc 48*59599516SKenneth E. Jansen nenl = lcblk(5,iblk) ! no. of vertices per element 49*59599516SKenneth E. Jansen iel = lcblk(1,iblk) 50*59599516SKenneth E. Jansen lelCat = lcblk(2,iblk) 51*59599516SKenneth E. Jansen lcsyst = lcblk(3,iblk) 52*59599516SKenneth E. Jansen iorder = lcblk(4,iblk) 53*59599516SKenneth E. Jansen nenl = lcblk(5,iblk) ! no. of vertices per element 54*59599516SKenneth E. Jansen nshl = lcblk(10,iblk) 55*59599516SKenneth E. Jansen mattyp = lcblk(7,iblk) 56*59599516SKenneth E. Jansen ndofl = lcblk(8,iblk) 57*59599516SKenneth E. Jansen nsymdl = lcblk(9,iblk) 58*59599516SKenneth E. Jansen npro = lcblk(1,iblk+1) - iel 59*59599516SKenneth E. Jansen ngauss = nint(lcsyst) 60*59599516SKenneth E. Jansenc 61*59599516SKenneth E. Jansenc.... compute and assemble the constarint factor, and mass matrix 62*59599516SKenneth E. Jansenc 63*59599516SKenneth E. Jansen 64*59599516SKenneth E. Jansen allocate (tmpshp(nshl,MAXQPT)) 65*59599516SKenneth E. Jansen allocate (tmpshgl(nsd,nshl,MAXQPT)) 66*59599516SKenneth E. Jansen 67*59599516SKenneth E. Jansen tmpshp(1:nshl,:) = shp(lcsyst,1:nshl,:) 68*59599516SKenneth E. Jansen tmpshgl(:,1:nshl,:) = shgl(lcsyst,:,1:nshl,:) 69*59599516SKenneth E. Jansenc 70*59599516SKenneth E. Jansen call volcon (y, x, tmpshp, 71*59599516SKenneth E. Jansen & tmpshgl, mien(iblk)%p, rmass, 72*59599516SKenneth E. Jansen & v_lambda1, hprime, v_lambda2) 73*59599516SKenneth E. Jansen 74*59599516SKenneth E. Jansen deallocate ( tmpshp ) 75*59599516SKenneth E. Jansen deallocate ( tmpshgl ) 76*59599516SKenneth E. Jansen enddo 77*59599516SKenneth E. Jansen 78*59599516SKenneth E. Jansenc 79*59599516SKenneth E. Jansenc ... multiple processor communication 80*59599516SKenneth E. Jansenc 81*59599516SKenneth E. Jansen if (numpe > 1) then 82*59599516SKenneth E. Jansen call commu (v_lambda1 , ilwork, 1 , 'in ') 83*59599516SKenneth E. Jansen call commu (v_lambda2 , ilwork, 1 , 'in ') 84*59599516SKenneth E. Jansen call commu (hprime , ilwork, 1 , 'in ') 85*59599516SKenneth E. Jansen call commu (rmass , ilwork, 1 , 'in ') 86*59599516SKenneth E. Jansen endif 87*59599516SKenneth E. Jansenc 88*59599516SKenneth E. Jansenc.... take care of periodic boundary conditions 89*59599516SKenneth E. Jansenc 90*59599516SKenneth E. Jansen do j= 1,nshg 91*59599516SKenneth E. Jansen if (btest(iBC(j),10)) then 92*59599516SKenneth E. Jansen i = iper(j) 93*59599516SKenneth E. Jansen rmass(i) = rmass(i) + rmass(j) 94*59599516SKenneth E. Jansen v_lambda1(i) = v_lambda1(i) + v_lambda1(j) 95*59599516SKenneth E. Jansen v_lambda2(i) = v_lambda2(i) + v_lambda2(j) 96*59599516SKenneth E. Jansen hprime(i) = hprime(i) + hprime(j) 97*59599516SKenneth E. Jansen endif 98*59599516SKenneth E. Jansen enddo 99*59599516SKenneth E. Jansenc 100*59599516SKenneth E. Jansen do j= 1,nshg 101*59599516SKenneth E. Jansen if (btest(iBC(j),10)) then 102*59599516SKenneth E. Jansen i = iper(j) 103*59599516SKenneth E. Jansen rmass(j) = rmass(i) 104*59599516SKenneth E. Jansen v_lambda1(j) = v_lambda1(i) 105*59599516SKenneth E. Jansen v_lambda2(j) = v_lambda2(i) 106*59599516SKenneth E. Jansen hprime(j) = hprime(i) 107*59599516SKenneth E. Jansen endif 108*59599516SKenneth E. Jansen enddo 109*59599516SKenneth E. Jansenc 110*59599516SKenneth E. Jansenc ... calculation of constraint factor 111*59599516SKenneth E. Jansenc 112*59599516SKenneth E. Jansen rmass = one/rmass 113*59599516SKenneth E. Jansen v_lambda1 = v_lambda1*rmass ! numerator of lambda 114*59599516SKenneth E. Jansen v_lambda2 = v_lambda2*rmass ! denominator of lambda 115*59599516SKenneth E. Jansen v_lambda = v_lambda1/(v_lambda2+epsM**2) 116*59599516SKenneth E. Jansen hprime = hprime*rmass 117*59599516SKenneth E. Jansen v_lambda = v_lambda*hprime 118*59599516SKenneth E. Jansenc 119*59599516SKenneth E. Jansenc ... commu out for the multiple processor 120*59599516SKenneth E. Jansenc 121*59599516SKenneth E. Jansen if(numpe > 1) then 122*59599516SKenneth E. Jansen call commu (v_lambda, ilwork, 1, 'out') 123*59599516SKenneth E. Jansen endif 124*59599516SKenneth E. Jansenc 125*59599516SKenneth E. Jansenc ... the following commented lines are for the different way of getting 126*59599516SKenneth E. Jansenc the denominator of constraint (lambda) calculation 127*59599516SKenneth E. Jansenc$$$ hprime=zero 128*59599516SKenneth E. Jansenc$$$ do kk=1, nshg 129*59599516SKenneth E. Jansenc$$$ if (abs (y(kk,6)) .le. epsilon_ls) then 130*59599516SKenneth E. Jansenc$$$ hprime(kk) = (0.5/epsilon_ls) * (1 131*59599516SKenneth E. Jansenc$$$ & + cos(pi*y(kk,6)/epsilon_ls)) 132*59599516SKenneth E. Jansenc$$$ endif 133*59599516SKenneth E. Jansenc$$$ enddo 134*59599516SKenneth E. Jansenc$$$ y(:,7) = y(:,7)+v_lambda*hprime/dtgl 135*59599516SKenneth E. Jansenc 136*59599516SKenneth E. Jansenc ... the vlome constraint applied on the second scalar 137*59599516SKenneth E. Jansenc 138*59599516SKenneth E. Jansen y(:,7) = y(:,7)+v_lambda/dtgl 139*59599516SKenneth E. Jansenc 140*59599516SKenneth E. Jansen return 141*59599516SKenneth E. Jansen end 142*59599516SKenneth E. Jansenc 143*59599516SKenneth E. Jansenc 144*59599516SKenneth E. Jansenc 145*59599516SKenneth E. Jansen 146*59599516SKenneth E. Jansen subroutine volcon (y, x, shp, 147*59599516SKenneth E. Jansen & shgl, ien, rmass, 148*59599516SKenneth E. Jansen & v_lambda1, hprime, v_lambda2) 149*59599516SKenneth E. Jansen 150*59599516SKenneth E. Jansenc--------------------------------------------------------------------- 151*59599516SKenneth E. Jansenc 152*59599516SKenneth E. Jansenc This subroutine is to calculate the element contribution to the 153*59599516SKenneth E. Jansenc constraint factor and mass matrix. 154*59599516SKenneth E. Jansenc 155*59599516SKenneth E. Jansenc--------------------------------------------------------------------- 156*59599516SKenneth E. Jansen include "common.h" 157*59599516SKenneth E. Jansenc 158*59599516SKenneth E. Jansen dimension y(nshg,ndof), x(numnp,nsd), 159*59599516SKenneth E. Jansen & shp(nshl,maxsh), 160*59599516SKenneth E. Jansen & shgl(nsd,nshl,maxsh), 161*59599516SKenneth E. Jansen & ien(npro,nshl), 162*59599516SKenneth E. Jansen & qres(nshg,idflx), rmass(nshg) 163*59599516SKenneth E. Jansenc 164*59599516SKenneth E. Jansenc.... element level declarations 165*59599516SKenneth E. Jansenc 166*59599516SKenneth E. Jansen dimension ycl(npro,nshl,ndof), xl(npro,nenl,nsd), 167*59599516SKenneth E. Jansen & rmassl(npro,nshl) 168*59599516SKenneth E. Jansen dimension sgn(npro,nshape), v_lambdal1(npro,nshl), 169*59599516SKenneth E. Jansen & v_lambda1(nshg), hprimel(npro,nshl), 170*59599516SKenneth E. Jansen & hprime(nshg), v_lambdal2(npro,nshl), 171*59599516SKenneth E. Jansen & v_lambda2(nshg) 172*59599516SKenneth E. Jansenc 173*59599516SKenneth E. Jansenc local arrays 174*59599516SKenneth E. Jansenc 175*59599516SKenneth E. Jansen dimension shg(npro,nshl,nsd), 176*59599516SKenneth E. Jansen & dxidx(npro,nsd,nsd), WdetJ(npro) 177*59599516SKenneth E. Jansenc 178*59599516SKenneth E. Jansen dimension shape(npro,nshl), 179*59599516SKenneth E. Jansen & shdrv(npro,nsd,nshl) 180*59599516SKenneth E. Jansenc 181*59599516SKenneth E. Jansenc 182*59599516SKenneth E. Jansenc.... for volume constraint calculation of redistancing step 183*59599516SKenneth E. Jansenc 184*59599516SKenneth E. Jansen dimension Sclr(npro), Sclrtmp(npro), 185*59599516SKenneth E. Jansen & h_prime(npro), tmp1(npro), 186*59599516SKenneth E. Jansen & tmp2(npro), 187*59599516SKenneth E. Jansen & v_lambdatmp(npro), v_lambdal(npro,nshl) 188*59599516SKenneth E. Jansenc$$$ & ,hprimel(npro,nshl), v_lambdal1(npro,nshl), 189*59599516SKenneth E. Jansenc$$$ & v_lambdal2(npro,nshl) 190*59599516SKenneth E. Jansenc above arrays must be uncommented for alternate method included below (commented) 191*59599516SKenneth E. Jansen real epsilon_tmp 192*59599516SKenneth E. Jansen 193*59599516SKenneth E. Jansenc 194*59599516SKenneth E. Jansenc.... create the matrix of mode signs for the hierarchic basis 195*59599516SKenneth E. Jansenc functions. 196*59599516SKenneth E. Jansenc 197*59599516SKenneth E. Jansen if (ipord .gt. 1) then 198*59599516SKenneth E. Jansen call getsgn(ien,sgn) 199*59599516SKenneth E. Jansen endif 200*59599516SKenneth E. Jansenc 201*59599516SKenneth E. Jansenc.... gather the variables 202*59599516SKenneth E. Jansenc 203*59599516SKenneth E. Jansen 204*59599516SKenneth E. Jansen call localy(y, ycl, ien, ndof, 'gather ') 205*59599516SKenneth E. Jansen call localx(x, xl, ien, nsd, 'gather ') 206*59599516SKenneth E. Jansenc 207*59599516SKenneth E. Jansenc.... get the element contributions of the numerator and denominator 208*59599516SKenneth E. Jansenc of lambda 209*59599516SKenneth E. Jansenc 210*59599516SKenneth E. Jansen rmassl = zero 211*59599516SKenneth E. Jansen v_lambdal1= zero 212*59599516SKenneth E. Jansen v_lambdal2= zero 213*59599516SKenneth E. Jansen hprimel = zero 214*59599516SKenneth E. Jansenc 215*59599516SKenneth E. Jansenc.... loop through the integration points 216*59599516SKenneth E. Jansenc 217*59599516SKenneth E. Jansen do intp = 1, ngauss 218*59599516SKenneth E. Jansen if (Qwt(lcsyst,intp) .eq. zero) cycle ! precaution 219*59599516SKenneth E. Jansenc 220*59599516SKenneth E. Jansenc.... create a matrix of shape functions (and derivatives) for each 221*59599516SKenneth E. Jansenc element at this quadrature point. These arrays will contain 222*59599516SKenneth E. Jansenc the correct signs for the hierarchic basis 223*59599516SKenneth E. Jansenc 224*59599516SKenneth E. Jansen call getshp(shp, shgl, sgn, 225*59599516SKenneth E. Jansen & shape, shdrv) 226*59599516SKenneth E. Jansenc 227*59599516SKenneth E. Jansenc.... initialize 228*59599516SKenneth E. Jansenc 229*59599516SKenneth E. Jansen h_prime = zero 230*59599516SKenneth E. Jansen sclr = zero 231*59599516SKenneth E. Jansen sclrtmp = zero 232*59599516SKenneth E. Jansen v_lambdatmp=zero 233*59599516SKenneth E. Jansen tmp1 =zero 234*59599516SKenneth E. Jansen tmp2 =zero 235*59599516SKenneth E. Jansen 236*59599516SKenneth E. Jansenc 237*59599516SKenneth E. Jansenc.... ---------------------> Element Metrics <----------------------- 238*59599516SKenneth E. Jansenc 239*59599516SKenneth E. Jansen call e3metric( xl, shdrv, dxidx, 240*59599516SKenneth E. Jansen & shg, WdetJ) 241*59599516SKenneth E. Jansen 242*59599516SKenneth E. Jansenc 243*59599516SKenneth E. Jansen do i = 1, nshl 244*59599516SKenneth E. Jansenc 245*59599516SKenneth E. Jansenc y(intp)=SUM_{a=1}^nshl (N_a(intp) Ya) 246*59599516SKenneth E. Jansenc 247*59599516SKenneth E. Jansen Sclr = Sclr + shape(:,i) * ycl(:,i,7) !d^kbar 248*59599516SKenneth E. Jansen sclrtmp = sclrtmp + shape(:,i) * ycl(:,i,6) !d^0 249*59599516SKenneth E. Jansen enddo 250*59599516SKenneth E. Jansen 251*59599516SKenneth E. Jansen if (isclr .eq. 2) then 252*59599516SKenneth E. Jansen epsilon_tmp = epsilon_lsd 253*59599516SKenneth E. Jansen else 254*59599516SKenneth E. Jansen epsilon_tmp = epsilon_ls 255*59599516SKenneth E. Jansen endif 256*59599516SKenneth E. Jansen 257*59599516SKenneth E. Jansen do i=1,npro 258*59599516SKenneth E. Jansen if (abs (Sclrtmp(i)) .le. epsilon_tmp) then 259*59599516SKenneth E. Jansen h_prime(i) = (0.5/epsilon_tmp) * (1 260*59599516SKenneth E. Jansen & + cos(pi*Sclrtmp(i)/epsilon_tmp)) 261*59599516SKenneth E. Jansen tmp1(i)=-h_prime(i)*(sclr(i)-sclrtmp(i))*dtgl 262*59599516SKenneth E. Jansen tmp2(i)=h_prime(i)**2 263*59599516SKenneth E. Jansenc v_lambdatmp(i)=tmp1(i)/(tmp2(i)+ epsM) 264*59599516SKenneth E. Jansen endif 265*59599516SKenneth E. Jansen enddo 266*59599516SKenneth E. Jansenc 267*59599516SKenneth E. Jansen do i=1,nshl 268*59599516SKenneth E. Jansen v_lambdal1(:,i)=v_lambdal1(:,i)+ shape(:,i)*WdetJ*tmp1 269*59599516SKenneth E. Jansen v_lambdal2(:,i)=v_lambdal2(:,i)+ shape(:,i)*WdetJ*tmp2 270*59599516SKenneth E. Jansenc v_lambdal(:,i)=v_lambdal(:,i)+ shape(:,i)*WdetJ*v_lambdatmp 271*59599516SKenneth E. Jansen hprimel(:,i) = hprimel(:,i) + shape(:,i)*WdetJ*h_prime 272*59599516SKenneth E. Jansen rmassl(:,i) =rmassl(:,i) + shape(:,i)*WdetJ 273*59599516SKenneth E. Jansen enddo 274*59599516SKenneth E. Jansenc 275*59599516SKenneth E. Jansenc.... end of the loop over integration points 276*59599516SKenneth E. Jansenc 277*59599516SKenneth E. Jansen enddo 278*59599516SKenneth E. Jansenc 279*59599516SKenneth E. Jansen call local (v_lambda1, v_lambdal1, ien, 1, 'scatter ') 280*59599516SKenneth E. Jansen call local (v_lambda2, v_lambdal2, ien, 1, 'scatter ') 281*59599516SKenneth E. Jansen call local (hprime, hprimel, ien, 1, 'scatter ') 282*59599516SKenneth E. Jansen call local (rmass, rmassl, ien, 1, 'scatter ') 283*59599516SKenneth E. Jansenc 284*59599516SKenneth E. Jansen return 285*59599516SKenneth E. Jansen end 286*59599516SKenneth E. Jansen 287