| 1 | c-----------------------------------------------------------------------
|
|---|
| 2 | subroutine gen_fast(df,sr,ss,st,x,y,z)
|
|---|
| 3 | c
|
|---|
| 4 | c Generate fast diagonalization matrices for each element
|
|---|
| 5 | c
|
|---|
| 6 | include 'SIZE'
|
|---|
| 7 | include 'INPUT'
|
|---|
| 8 | include 'PARALLEL'
|
|---|
| 9 | include 'SOLN'
|
|---|
| 10 | include 'WZ'
|
|---|
| 11 | c
|
|---|
| 12 | parameter(lxx=lx1*lx1)
|
|---|
| 13 | real df(lx1*ly1*lz1,1),sr(lxx*2,1),ss(lxx*2,1),st(lxx*2,1)
|
|---|
| 14 | c
|
|---|
| 15 | common /ctmpf/ lr(2*lx1+4),ls(2*lx1+4),lt(2*lx1+4)
|
|---|
| 16 | $ , llr(lelt),lls(lelt),llt(lelt)
|
|---|
| 17 | $ , lmr(lelt),lms(lelt),lmt(lelt)
|
|---|
| 18 | $ , lrr(lelt),lrs(lelt),lrt(lelt)
|
|---|
| 19 | real lr ,ls ,lt
|
|---|
| 20 | real llr,lls,llt
|
|---|
| 21 | real lmr,lms,lmt
|
|---|
| 22 | real lrr,lrs,lrt
|
|---|
| 23 | c
|
|---|
| 24 | integer lbr,rbr,lbs,rbs,lbt,rbt,e
|
|---|
| 25 | c
|
|---|
| 26 | real x(lx1,ly1,lz1,nelv)
|
|---|
| 27 | real y(lx1,ly1,lz1,nelv)
|
|---|
| 28 | real z(lx1,ly1,lz1,nelv)
|
|---|
| 29 | real axwt(lx2)
|
|---|
| 30 |
|
|---|
| 31 | ierr = 0
|
|---|
| 32 |
|
|---|
| 33 | do e=1,nelv
|
|---|
| 34 | c
|
|---|
| 35 | if (param(44).eq.1) then
|
|---|
| 36 | call get_fast_bc(lbr,rbr,lbs,rbs,lbt,rbt,e,2,ierr)
|
|---|
| 37 | else
|
|---|
| 38 | call get_fast_bc(lbr,rbr,lbs,rbs,lbt,rbt,e,3,ierr)
|
|---|
| 39 | endif
|
|---|
| 40 | c
|
|---|
| 41 | c Set up matrices for each element.
|
|---|
| 42 | c
|
|---|
| 43 | if (param(44).eq.1) then
|
|---|
| 44 | call set_up_fast_1D_fem( sr(1,e),lr,nr ,lbr,rbr
|
|---|
| 45 | $ ,llr(e),lmr(e),lrr(e),zgm2(1,1),lx2,e)
|
|---|
| 46 | else
|
|---|
| 47 | call set_up_fast_1D_sem( sr(1,e),lr,nr ,lbr,rbr
|
|---|
| 48 | $ ,llr(e),lmr(e),lrr(e),e)
|
|---|
| 49 | endif
|
|---|
| 50 | if (ifaxis) then
|
|---|
| 51 | xsum = vlsum(wxm2,lx2)
|
|---|
| 52 | do i=1,ly2
|
|---|
| 53 | yavg = vlsc2(y(1,i,1,e),wxm2,lx2)/xsum
|
|---|
| 54 | axwt(i) = yavg
|
|---|
| 55 | enddo
|
|---|
| 56 | call set_up_fast_1D_fem_ax( ss(1,e),ls,ns ,lbs,rbs
|
|---|
| 57 | $ ,lls(e),lms(e),lrs(e),zgm2(1,2),axwt,ly2,e)
|
|---|
| 58 | else
|
|---|
| 59 | if (param(44).eq.1) then
|
|---|
| 60 | call set_up_fast_1D_fem( ss(1,e),ls,ns ,lbs,rbs
|
|---|
| 61 | $ ,lls(e),lms(e),lrs(e),zgm2(1,2),ly2,e)
|
|---|
| 62 | else
|
|---|
| 63 | call set_up_fast_1D_sem( ss(1,e),ls,ns ,lbs,rbs
|
|---|
| 64 | $ ,lls(e),lms(e),lrs(e),e)
|
|---|
| 65 | endif
|
|---|
| 66 | endif
|
|---|
| 67 | if (if3d) then
|
|---|
| 68 | if (param(44).eq.1) then
|
|---|
| 69 | call set_up_fast_1D_fem( st(1,e),lt,nt ,lbt,rbt
|
|---|
| 70 | $ ,llt(e),lmt(e),lrt(e),zgm2(1,3),lz2,e)
|
|---|
| 71 | else
|
|---|
| 72 | call set_up_fast_1D_sem( st(1,e),lt,nt ,lbt,rbt
|
|---|
| 73 | $ ,llt(e),lmt(e),lrt(e),e)
|
|---|
| 74 | endif
|
|---|
| 75 | endif
|
|---|
| 76 | c
|
|---|
| 77 | c DIAGNOSTICS
|
|---|
| 78 | c
|
|---|
| 79 | c n12 = min(9,nr)
|
|---|
| 80 | c write(6,1) e,'1D lr',llr(e),lmr(e),lrr(e),(lr(k),k=1,n12)
|
|---|
| 81 | c write(6,1) e,'1D ls',lls(e),lms(e),lrs(e),(ls(k),k=1,n12)
|
|---|
| 82 | c if (if3d)
|
|---|
| 83 | c $ write(6,1) e,'1D lt',llt(e),lmt(e),lrt(e),(lt(k),k=1,n12)
|
|---|
| 84 | c 1 format(i6,1x,a5,1p12e12.4)
|
|---|
| 85 | c
|
|---|
| 86 | c
|
|---|
| 87 | c Set up diagonal inverse
|
|---|
| 88 | c
|
|---|
| 89 | if (if3d) then
|
|---|
| 90 | eps = 1.e-5 * (vlmax(lr(2),nr-2)
|
|---|
| 91 | $ + vlmax(ls(2),ns-2) + vlmax(lt(2),nt-2))
|
|---|
| 92 | l = 1
|
|---|
| 93 | do k=1,nt
|
|---|
| 94 | do j=1,ns
|
|---|
| 95 | do i=1,nr
|
|---|
| 96 | diag = lr(i) + ls(j) + lt(k)
|
|---|
| 97 | if (diag.gt.eps) then
|
|---|
| 98 | df(l,e) = 1.0/diag
|
|---|
| 99 | else
|
|---|
| 100 | c write(6,3) e,'Reset Eig in gen fast:',i,j,k,l
|
|---|
| 101 | c $ ,eps,diag,lr(i),ls(j),lt(k)
|
|---|
| 102 | c 3 format(i6,1x,a21,4i5,1p5e12.4)
|
|---|
| 103 | df(l,e) = 0.0
|
|---|
| 104 | endif
|
|---|
| 105 | l = l+1
|
|---|
| 106 | enddo
|
|---|
| 107 | enddo
|
|---|
| 108 | enddo
|
|---|
| 109 | else
|
|---|
| 110 | eps = 1.e-5*(vlmax(lr(2),nr-2) + vlmax(ls(2),ns-2))
|
|---|
| 111 | l = 1
|
|---|
| 112 | do j=1,ns
|
|---|
| 113 | do i=1,nr
|
|---|
| 114 | diag = lr(i) + ls(j)
|
|---|
| 115 | if (diag.gt.eps) then
|
|---|
| 116 | df(l,e) = 1.0/diag
|
|---|
| 117 | else
|
|---|
| 118 | c write(6,2) e,'Reset Eig in gen fast:',i,j,l
|
|---|
| 119 | c $ ,eps,diag,lr(i),ls(j)
|
|---|
| 120 | c 2 format(i6,1x,a21,3i5,1p4e12.4)
|
|---|
| 121 | df(l,e) = 0.0
|
|---|
| 122 | endif
|
|---|
| 123 | l = l+1
|
|---|
| 124 | enddo
|
|---|
| 125 | enddo
|
|---|
| 126 | endif
|
|---|
| 127 | c
|
|---|
| 128 | c Next element ....
|
|---|
| 129 | c
|
|---|
| 130 | enddo
|
|---|
| 131 |
|
|---|
| 132 | ierrmx = iglmax(ierr,1)
|
|---|
| 133 | if (ierrmx.gt.0) then
|
|---|
| 134 | if (ierr.gt.0) write(6,*) nid,ierr,' BC FAIL'
|
|---|
| 135 | call exitti('E INVALID BC FOUND in genfast$',ierrmx)
|
|---|
| 136 | endif
|
|---|
| 137 |
|
|---|
| 138 |
|
|---|
| 139 | return
|
|---|
| 140 | end
|
|---|
| 141 | c-----------------------------------------------------------------------
|
|---|
| 142 | subroutine gen_fast_spacing(x,y,z)
|
|---|
| 143 | c
|
|---|
| 144 | c Generate fast diagonalization matrices for each element
|
|---|
| 145 | c
|
|---|
| 146 | include 'SIZE'
|
|---|
| 147 | include 'INPUT'
|
|---|
| 148 | include 'PARALLEL'
|
|---|
| 149 | include 'SOLN'
|
|---|
| 150 | include 'WZ'
|
|---|
| 151 | c
|
|---|
| 152 | parameter(lxx=lx1*lx1)
|
|---|
| 153 | c
|
|---|
| 154 | common /ctmpf/ lr(2*lx1+4),ls(2*lx1+4),lt(2*lx1+4)
|
|---|
| 155 | $ , llr(lelt),lls(lelt),llt(lelt)
|
|---|
| 156 | $ , lmr(lelt),lms(lelt),lmt(lelt)
|
|---|
| 157 | $ , lrr(lelt),lrs(lelt),lrt(lelt)
|
|---|
| 158 | real lr ,ls ,lt
|
|---|
| 159 | real llr,lls,llt
|
|---|
| 160 | real lmr,lms,lmt
|
|---|
| 161 | real lrr,lrs,lrt
|
|---|
| 162 | c
|
|---|
| 163 | integer lbr,rbr,lbs,rbs,lbt,rbt,e
|
|---|
| 164 | c
|
|---|
| 165 | real x(lx1,ly1,lz1,nelv)
|
|---|
| 166 | real y(lx1,ly1,lz1,nelv)
|
|---|
| 167 | real z(lx1,ly1,lz1,nelv)
|
|---|
| 168 | real axwt(lx2)
|
|---|
| 169 |
|
|---|
| 170 | ierr = 0
|
|---|
| 171 |
|
|---|
| 172 | if (param(44).eq.1) then
|
|---|
| 173 | c __ __ __
|
|---|
| 174 | c Now, for each element, compute lr,ls,lt between specified planes
|
|---|
| 175 | c
|
|---|
| 176 | n1 = lx2
|
|---|
| 177 | n2 = lx2+1
|
|---|
| 178 | nz0 = 1
|
|---|
| 179 | nzn = 1
|
|---|
| 180 | if (if3d) then
|
|---|
| 181 | nz0= 0
|
|---|
| 182 | nzn=n2
|
|---|
| 183 | endif
|
|---|
| 184 | eps = 1.e-7
|
|---|
| 185 | if (wdsize.eq.8) eps = 1.e-14
|
|---|
| 186 | c
|
|---|
| 187 | c Find mean spacing between "left-most" planes
|
|---|
| 188 | call plane_space2(llr,lls,llt, 0,wxm2,x,y,z,n1,n2,nz0,nzn)
|
|---|
| 189 | c
|
|---|
| 190 | c Find mean spacing between "middle" planes
|
|---|
| 191 | call plane_space (lmr,lms,lmt, 1,n1,wxm2,x,y,z,n1,n2,nz0,nzn)
|
|---|
| 192 | c
|
|---|
| 193 | c Find mean spacing between "right-most" planes
|
|---|
| 194 | call plane_space2(lrr,lrs,lrt,n2,wxm2,x,y,z,n1,n2,nz0,nzn)
|
|---|
| 195 | c
|
|---|
| 196 | else
|
|---|
| 197 | call load_semhat_weighted ! Fills the SEMHAT arrays
|
|---|
| 198 | endif
|
|---|
| 199 |
|
|---|
| 200 | return
|
|---|
| 201 | end
|
|---|
| 202 | c-----------------------------------------------------------------------
|
|---|
| 203 | subroutine plane_space_std(lr,ls,lt,i1,i2,w,x,y,z,nx,nxn,nz0,nzn)
|
|---|
| 204 | c
|
|---|
| 205 | c This routine now replaced by "plane_space()"
|
|---|
| 206 | c
|
|---|
| 207 | c Here, spacing is based on arithmetic mean.
|
|---|
| 208 | c New verision uses harmonic mean. pff 2/10/07
|
|---|
| 209 | c
|
|---|
| 210 | include 'SIZE'
|
|---|
| 211 | include 'INPUT'
|
|---|
| 212 | c
|
|---|
| 213 | real w(1),lr(1),ls(1),lt(1)
|
|---|
| 214 | real x(0:nxn,0:nxn,nz0:nzn,1)
|
|---|
| 215 | real y(0:nxn,0:nxn,nz0:nzn,1)
|
|---|
| 216 | real z(0:nxn,0:nxn,nz0:nzn,1)
|
|---|
| 217 | real lr2,ls2,lt2
|
|---|
| 218 | c __ __ __
|
|---|
| 219 | c Now, for each element, compute lr,ls,lt between specified planes
|
|---|
| 220 | c
|
|---|
| 221 | ny = nx
|
|---|
| 222 | nz = nx
|
|---|
| 223 | j1 = i1
|
|---|
| 224 | k1 = i1
|
|---|
| 225 | j2 = i2
|
|---|
| 226 | k2 = i2
|
|---|
| 227 | c
|
|---|
| 228 | do ie=1,nelv
|
|---|
| 229 | c
|
|---|
| 230 | if (if3d) then
|
|---|
| 231 | lr2 = 0.
|
|---|
| 232 | wsum = 0.
|
|---|
| 233 | do k=1,nz
|
|---|
| 234 | do j=1,ny
|
|---|
| 235 | weight = w(j)*w(k)
|
|---|
| 236 | lr2 = lr2 + ( (x(i2,j,k,ie)-x(i1,j,k,ie))**2
|
|---|
| 237 | $ + (y(i2,j,k,ie)-y(i1,j,k,ie))**2
|
|---|
| 238 | $ + (z(i2,j,k,ie)-z(i1,j,k,ie))**2 )
|
|---|
| 239 | $ * weight
|
|---|
| 240 | wsum = wsum + weight
|
|---|
| 241 | enddo
|
|---|
| 242 | enddo
|
|---|
| 243 | lr2 = lr2/wsum
|
|---|
| 244 | lr(ie) = sqrt(lr2)
|
|---|
| 245 | c
|
|---|
| 246 | ls2 = 0.
|
|---|
| 247 | wsum = 0.
|
|---|
| 248 | do k=1,nz
|
|---|
| 249 | do i=1,nx
|
|---|
| 250 | weight = w(i)*w(k)
|
|---|
| 251 | ls2 = ls2 + ( (x(i,j2,k,ie)-x(i,j1,k,ie))**2
|
|---|
| 252 | $ + (y(i,j2,k,ie)-y(i,j1,k,ie))**2
|
|---|
| 253 | $ + (z(i,j2,k,ie)-z(i,j1,k,ie))**2 )
|
|---|
| 254 | $ * weight
|
|---|
| 255 | wsum = wsum + weight
|
|---|
| 256 | enddo
|
|---|
| 257 | enddo
|
|---|
| 258 | ls2 = ls2/wsum
|
|---|
| 259 | ls(ie) = sqrt(ls2)
|
|---|
| 260 | c
|
|---|
| 261 | lt2 = 0.
|
|---|
| 262 | wsum = 0.
|
|---|
| 263 | do j=1,ny
|
|---|
| 264 | do i=1,nx
|
|---|
| 265 | weight = w(i)*w(j)
|
|---|
| 266 | lt2 = lt2 + ( (x(i,j,k2,ie)-x(i,j,k1,ie))**2
|
|---|
| 267 | $ + (y(i,j,k2,ie)-y(i,j,k1,ie))**2
|
|---|
| 268 | $ + (z(i,j,k2,ie)-z(i,j,k1,ie))**2 )
|
|---|
| 269 | $ * weight
|
|---|
| 270 | wsum = wsum + weight
|
|---|
| 271 | enddo
|
|---|
| 272 | enddo
|
|---|
| 273 | lt2 = lt2/wsum
|
|---|
| 274 | lt(ie) = sqrt(lt2)
|
|---|
| 275 | c
|
|---|
| 276 | else
|
|---|
| 277 | lr2 = 0.
|
|---|
| 278 | wsum = 0.
|
|---|
| 279 | do j=1,ny
|
|---|
| 280 | weight = w(j)
|
|---|
| 281 | lr2 = lr2 + ( (x(i2,j,1,ie)-x(i1,j,1,ie))**2
|
|---|
| 282 | $ + (y(i2,j,1,ie)-y(i1,j,1,ie))**2 )
|
|---|
| 283 | $ * weight
|
|---|
| 284 | wsum = wsum + weight
|
|---|
| 285 | enddo
|
|---|
| 286 | lr2 = lr2/wsum
|
|---|
| 287 | lr(ie) = sqrt(lr2)
|
|---|
| 288 | c
|
|---|
| 289 | ls2 = 0.
|
|---|
| 290 | wsum = 0.
|
|---|
| 291 | do i=1,nx
|
|---|
| 292 | weight = w(i)
|
|---|
| 293 | ls2 = ls2 + ( (x(i,j2,1,ie)-x(i,j1,1,ie))**2
|
|---|
| 294 | $ + (y(i,j2,1,ie)-y(i,j1,1,ie))**2 )
|
|---|
| 295 | $ * weight
|
|---|
| 296 | wsum = wsum + weight
|
|---|
| 297 | enddo
|
|---|
| 298 | ls2 = ls2/wsum
|
|---|
| 299 | ls(ie) = sqrt(ls2)
|
|---|
| 300 | c write(6,*) 'lrls',ie,lr(ie),ls(ie)
|
|---|
| 301 | endif
|
|---|
| 302 | enddo
|
|---|
| 303 | return
|
|---|
| 304 | end
|
|---|
| 305 | c-----------------------------------------------------------------------
|
|---|
| 306 | subroutine plane_space(lr,ls,lt,i1,i2,w,x,y,z,nx,nxn,nz0,nzn)
|
|---|
| 307 | c
|
|---|
| 308 | c Here, spacing is based on harmonic mean. pff 2/10/07
|
|---|
| 309 | c
|
|---|
| 310 | c
|
|---|
| 311 | include 'SIZE'
|
|---|
| 312 | include 'INPUT'
|
|---|
| 313 | c
|
|---|
| 314 | real w(1),lr(1),ls(1),lt(1)
|
|---|
| 315 | real x(0:nxn,0:nxn,nz0:nzn,1)
|
|---|
| 316 | real y(0:nxn,0:nxn,nz0:nzn,1)
|
|---|
| 317 | real z(0:nxn,0:nxn,nz0:nzn,1)
|
|---|
| 318 | real lr2,ls2,lt2
|
|---|
| 319 | c __ __ __
|
|---|
| 320 | c Now, for each element, compute lr,ls,lt between specified planes
|
|---|
| 321 | c
|
|---|
| 322 | ny = nx
|
|---|
| 323 | nz = nx
|
|---|
| 324 | j1 = i1
|
|---|
| 325 | k1 = i1
|
|---|
| 326 | j2 = i2
|
|---|
| 327 | k2 = i2
|
|---|
| 328 | c
|
|---|
| 329 | do ie=1,nelv
|
|---|
| 330 | c
|
|---|
| 331 | if (if3d) then
|
|---|
| 332 | lr2 = 0.
|
|---|
| 333 | wsum = 0.
|
|---|
| 334 | do k=1,nz
|
|---|
| 335 | do j=1,ny
|
|---|
| 336 | weight = w(j)*w(k)
|
|---|
| 337 | c lr2 = lr2 + ( (x(i2,j,k,ie)-x(i1,j,k,ie))**2
|
|---|
| 338 | c $ + (y(i2,j,k,ie)-y(i1,j,k,ie))**2
|
|---|
| 339 | c $ + (z(i2,j,k,ie)-z(i1,j,k,ie))**2 )
|
|---|
| 340 | c $ * weight
|
|---|
| 341 | lr2 = lr2 + weight /
|
|---|
| 342 | $ ( (x(i2,j,k,ie)-x(i1,j,k,ie))**2
|
|---|
| 343 | $ + (y(i2,j,k,ie)-y(i1,j,k,ie))**2
|
|---|
| 344 | $ + (z(i2,j,k,ie)-z(i1,j,k,ie))**2 )
|
|---|
| 345 | wsum = wsum + weight
|
|---|
| 346 | enddo
|
|---|
| 347 | enddo
|
|---|
| 348 | lr2 = lr2/wsum
|
|---|
| 349 | lr(ie) = 1./sqrt(lr2)
|
|---|
| 350 | c
|
|---|
| 351 | ls2 = 0.
|
|---|
| 352 | wsum = 0.
|
|---|
| 353 | do k=1,nz
|
|---|
| 354 | do i=1,nx
|
|---|
| 355 | weight = w(i)*w(k)
|
|---|
| 356 | c ls2 = ls2 + ( (x(i,j2,k,ie)-x(i,j1,k,ie))**2
|
|---|
| 357 | c $ + (y(i,j2,k,ie)-y(i,j1,k,ie))**2
|
|---|
| 358 | c $ + (z(i,j2,k,ie)-z(i,j1,k,ie))**2 )
|
|---|
| 359 | c $ * weight
|
|---|
| 360 | ls2 = ls2 + weight /
|
|---|
| 361 | $ ( (x(i,j2,k,ie)-x(i,j1,k,ie))**2
|
|---|
| 362 | $ + (y(i,j2,k,ie)-y(i,j1,k,ie))**2
|
|---|
| 363 | $ + (z(i,j2,k,ie)-z(i,j1,k,ie))**2 )
|
|---|
| 364 | wsum = wsum + weight
|
|---|
| 365 | enddo
|
|---|
| 366 | enddo
|
|---|
| 367 | ls2 = ls2/wsum
|
|---|
| 368 | ls(ie) = 1./sqrt(ls2)
|
|---|
| 369 | c
|
|---|
| 370 | lt2 = 0.
|
|---|
| 371 | wsum = 0.
|
|---|
| 372 | do j=1,ny
|
|---|
| 373 | do i=1,nx
|
|---|
| 374 | weight = w(i)*w(j)
|
|---|
| 375 | c lt2 = lt2 + ( (x(i,j,k2,ie)-x(i,j,k1,ie))**2
|
|---|
| 376 | c $ + (y(i,j,k2,ie)-y(i,j,k1,ie))**2
|
|---|
| 377 | c $ + (z(i,j,k2,ie)-z(i,j,k1,ie))**2 )
|
|---|
| 378 | c $ * weight
|
|---|
| 379 | lt2 = lt2 + weight /
|
|---|
| 380 | $ ( (x(i,j,k2,ie)-x(i,j,k1,ie))**2
|
|---|
| 381 | $ + (y(i,j,k2,ie)-y(i,j,k1,ie))**2
|
|---|
| 382 | $ + (z(i,j,k2,ie)-z(i,j,k1,ie))**2 )
|
|---|
| 383 | wsum = wsum + weight
|
|---|
| 384 | enddo
|
|---|
| 385 | enddo
|
|---|
| 386 | lt2 = lt2/wsum
|
|---|
| 387 | lt(ie) = 1./sqrt(lt2)
|
|---|
| 388 | c
|
|---|
| 389 | else ! 2D
|
|---|
| 390 | lr2 = 0.
|
|---|
| 391 | wsum = 0.
|
|---|
| 392 | do j=1,ny
|
|---|
| 393 | weight = w(j)
|
|---|
| 394 | c lr2 = lr2 + ( (x(i2,j,1,ie)-x(i1,j,1,ie))**2
|
|---|
| 395 | c $ + (y(i2,j,1,ie)-y(i1,j,1,ie))**2 )
|
|---|
| 396 | c $ * weight
|
|---|
| 397 | lr2 = lr2 + weight /
|
|---|
| 398 | $ ( (x(i2,j,1,ie)-x(i1,j,1,ie))**2
|
|---|
| 399 | $ + (y(i2,j,1,ie)-y(i1,j,1,ie))**2 )
|
|---|
| 400 | wsum = wsum + weight
|
|---|
| 401 | enddo
|
|---|
| 402 | lr2 = lr2/wsum
|
|---|
| 403 | lr(ie) = 1./sqrt(lr2)
|
|---|
| 404 | c
|
|---|
| 405 | ls2 = 0.
|
|---|
| 406 | wsum = 0.
|
|---|
| 407 | do i=1,nx
|
|---|
| 408 | weight = w(i)
|
|---|
| 409 | c ls2 = ls2 + ( (x(i,j2,1,ie)-x(i,j1,1,ie))**2
|
|---|
| 410 | c $ + (y(i,j2,1,ie)-y(i,j1,1,ie))**2 )
|
|---|
| 411 | c $ * weight
|
|---|
| 412 | ls2 = ls2 + weight /
|
|---|
| 413 | $ ( (x(i,j2,1,ie)-x(i,j1,1,ie))**2
|
|---|
| 414 | $ + (y(i,j2,1,ie)-y(i,j1,1,ie))**2 )
|
|---|
| 415 | wsum = wsum + weight
|
|---|
| 416 | enddo
|
|---|
| 417 | ls2 = ls2/wsum
|
|---|
| 418 | ls(ie) = 1./sqrt(ls2)
|
|---|
| 419 | c write(6,*) 'lrls',ie,lr(ie),ls(ie)
|
|---|
| 420 | endif
|
|---|
| 421 | enddo
|
|---|
| 422 | return
|
|---|
| 423 | end
|
|---|
| 424 | c-----------------------------------------------------------------------
|
|---|
| 425 | subroutine plane_space2(lr,ls,lt,i1,w,x,y,z,nx,nxn,nz0,nzn)
|
|---|
| 426 | c
|
|---|
| 427 | c Here, the local spacing is already given in the surface term.
|
|---|
| 428 | c This addition made to simplify the periodic bdry treatment.
|
|---|
| 429 | c
|
|---|
| 430 | c
|
|---|
| 431 | include 'SIZE'
|
|---|
| 432 | include 'INPUT'
|
|---|
| 433 | c
|
|---|
| 434 | real w(1),lr(1),ls(1),lt(1)
|
|---|
| 435 | real x(0:nxn,0:nxn,nz0:nzn,1)
|
|---|
| 436 | real y(0:nxn,0:nxn,nz0:nzn,1)
|
|---|
| 437 | real z(0:nxn,0:nxn,nz0:nzn,1)
|
|---|
| 438 | real lr2,ls2,lt2
|
|---|
| 439 | c __ __ __
|
|---|
| 440 | c Now, for each element, compute lr,ls,lt between specified planes
|
|---|
| 441 | c
|
|---|
| 442 | ny = nx
|
|---|
| 443 | nz = nx
|
|---|
| 444 | j1 = i1
|
|---|
| 445 | k1 = i1
|
|---|
| 446 | c
|
|---|
| 447 | do ie=1,nelv
|
|---|
| 448 | c
|
|---|
| 449 | if (if3d) then
|
|---|
| 450 | lr2 = 0.
|
|---|
| 451 | wsum = 0.
|
|---|
| 452 | do k=1,nz
|
|---|
| 453 | do j=1,ny
|
|---|
| 454 | weight = w(j)*w(k)
|
|---|
| 455 | lr2 = lr2 + ( (x(i1,j,k,ie))**2
|
|---|
| 456 | $ + (y(i1,j,k,ie))**2
|
|---|
| 457 | $ + (z(i1,j,k,ie))**2 )
|
|---|
| 458 | $ * weight
|
|---|
| 459 | wsum = wsum + weight
|
|---|
| 460 | enddo
|
|---|
| 461 | enddo
|
|---|
| 462 | lr2 = lr2/wsum
|
|---|
| 463 | lr(ie) = sqrt(lr2)
|
|---|
| 464 | c
|
|---|
| 465 | ls2 = 0.
|
|---|
| 466 | wsum = 0.
|
|---|
| 467 | do k=1,nz
|
|---|
| 468 | do i=1,nx
|
|---|
| 469 | weight = w(i)*w(k)
|
|---|
| 470 | ls2 = ls2 + ( (x(i,j1,k,ie))**2
|
|---|
| 471 | $ + (y(i,j1,k,ie))**2
|
|---|
| 472 | $ + (z(i,j1,k,ie))**2 )
|
|---|
| 473 | $ * weight
|
|---|
| 474 | wsum = wsum + weight
|
|---|
| 475 | enddo
|
|---|
| 476 | enddo
|
|---|
| 477 | ls2 = ls2/wsum
|
|---|
| 478 | ls(ie) = sqrt(ls2)
|
|---|
| 479 | c
|
|---|
| 480 | lt2 = 0.
|
|---|
| 481 | wsum = 0.
|
|---|
| 482 | do j=1,ny
|
|---|
| 483 | do i=1,nx
|
|---|
| 484 | weight = w(i)*w(j)
|
|---|
| 485 | lt2 = lt2 + ( (x(i,j,k1,ie))**2
|
|---|
| 486 | $ + (y(i,j,k1,ie))**2
|
|---|
| 487 | $ + (z(i,j,k1,ie))**2 )
|
|---|
| 488 | $ * weight
|
|---|
| 489 | wsum = wsum + weight
|
|---|
| 490 | enddo
|
|---|
| 491 | enddo
|
|---|
| 492 | lt2 = lt2/wsum
|
|---|
| 493 | lt(ie) = sqrt(lt2)
|
|---|
| 494 | c write(6,1) 'lrlslt',ie,lr(ie),ls(ie),lt(ie)
|
|---|
| 495 | 1 format(a6,i5,1p3e12.4)
|
|---|
| 496 | c
|
|---|
| 497 | else
|
|---|
| 498 | lr2 = 0.
|
|---|
| 499 | wsum = 0.
|
|---|
| 500 | do j=1,ny
|
|---|
| 501 | weight = w(j)
|
|---|
| 502 | lr2 = lr2 + ( (x(i1,j,1,ie))**2
|
|---|
| 503 | $ + (y(i1,j,1,ie))**2 )
|
|---|
| 504 | $ * weight
|
|---|
| 505 | wsum = wsum + weight
|
|---|
| 506 | enddo
|
|---|
| 507 | lr2 = lr2/wsum
|
|---|
| 508 | lr(ie) = sqrt(lr2)
|
|---|
| 509 | c
|
|---|
| 510 | ls2 = 0.
|
|---|
| 511 | wsum = 0.
|
|---|
| 512 | do i=1,nx
|
|---|
| 513 | weight = w(i)
|
|---|
| 514 | ls2 = ls2 + ( (x(i,j1,1,ie))**2
|
|---|
| 515 | $ + (y(i,j1,1,ie))**2 )
|
|---|
| 516 | $ * weight
|
|---|
| 517 | wsum = wsum + weight
|
|---|
| 518 | enddo
|
|---|
| 519 | ls2 = ls2/wsum
|
|---|
| 520 | ls(ie) = sqrt(ls2)
|
|---|
| 521 | c write(6,*) 'lrls',ie,lr(ie),ls(ie),lt(ie)
|
|---|
| 522 | endif
|
|---|
| 523 | enddo
|
|---|
| 524 | return
|
|---|
| 525 | end
|
|---|
| 526 | c-----------------------------------------------------------------------
|
|---|
| 527 | subroutine set_up_fast_1D_fem(s,lam,n,lbc,rbc,ll,lm,lr,z,nz,e)
|
|---|
| 528 | real s(1),lam(1),ll,lm,lr,z(1)
|
|---|
| 529 | integer lbc,rbc,e
|
|---|
| 530 | c
|
|---|
| 531 | parameter (m=100)
|
|---|
| 532 | real dx(0:m)
|
|---|
| 533 | integer icalld
|
|---|
| 534 | save icalld
|
|---|
| 535 | data icalld/0/
|
|---|
| 536 | c
|
|---|
| 537 | icalld=icalld+1
|
|---|
| 538 | c
|
|---|
| 539 | if (nz.gt.m-3) then
|
|---|
| 540 | write(6,*) 'ABORT. Error in set_up_fast_1D_fem. Increase m to'
|
|---|
| 541 | $ , nz
|
|---|
| 542 | call exitt
|
|---|
| 543 | endif
|
|---|
| 544 | c
|
|---|
| 545 | c In the present scheme, each element is viewed as a d-fold
|
|---|
| 546 | c tensor of (1+nz+1) arrays, even if funky bc's are applied
|
|---|
| 547 | c on either end of the 1D array.
|
|---|
| 548 | c
|
|---|
| 549 | n = nz+2
|
|---|
| 550 | c
|
|---|
| 551 | c Compute spacing, dx()
|
|---|
| 552 | c
|
|---|
| 553 | call set_up_1D_geom(dx,lbc,rbc,ll,lm,lr,z,nz)
|
|---|
| 554 | c
|
|---|
| 555 | nn1 = n*n + 1
|
|---|
| 556 | call gen_eigs_A_fem(s,s(nn1),lam,n,dx,lbc,rbc)
|
|---|
| 557 | c
|
|---|
| 558 | return
|
|---|
| 559 | end
|
|---|
| 560 | c-----------------------------------------------------------------------
|
|---|
| 561 | subroutine set_up_1D_geom(dx,lbc,rbc,ll,lm,lr,z,nz)
|
|---|
| 562 | c
|
|---|
| 563 | real dx(0:1),ll,lm,lr,z(2)
|
|---|
| 564 | integer lbc,rbc
|
|---|
| 565 | c
|
|---|
| 566 | c
|
|---|
| 567 | c Set up the 1D geometry for the tensor-product based overlapping Schwarz
|
|---|
| 568 | c
|
|---|
| 569 | c Upon return:
|
|---|
| 570 | c
|
|---|
| 571 | c dx() contains the spacing required to set up the stiffness matrix.
|
|---|
| 572 | c
|
|---|
| 573 | c
|
|---|
| 574 | c Input:
|
|---|
| 575 | c
|
|---|
| 576 | c lbc (rbc) is 0 if the left (right) BC is Dirichlet, 1 if Neumann.
|
|---|
| 577 | c
|
|---|
| 578 | c ll is the space between the left-most Gauss point of the middle
|
|---|
| 579 | c element and the right-most Gauss point of the LEFT element
|
|---|
| 580 | c
|
|---|
| 581 | c lm is the space between the left-most Gauss point of the middle
|
|---|
| 582 | c element and the right-most Gauss point of the MIDDLE element
|
|---|
| 583 | c
|
|---|
| 584 | c lr is the space between the right-most Gauss point of the middle
|
|---|
| 585 | c element and the left-most Gauss point of the RIGHT element
|
|---|
| 586 | c
|
|---|
| 587 | c --- if ll (lr) is very small (0), it indicates that there is no
|
|---|
| 588 | c left (right) spacing, and that the left (right) BC is Neumann.
|
|---|
| 589 | c
|
|---|
| 590 | c
|
|---|
| 591 | c z() is the array of nz Gauss points on the interval ]-1,1[.
|
|---|
| 592 | c
|
|---|
| 593 | c Boundary conditions:
|
|---|
| 594 | c
|
|---|
| 595 | c bc = 0 -- std. Dirichlet bc applied 2 points away from interior
|
|---|
| 596 | c bc = 1 -- Dirichlet bc applied 1 point away from interior (outflow)
|
|---|
| 597 | c bc = 2 -- Neumann bc applied on interior point (W,v,V,SYM,...)
|
|---|
| 598 | c
|
|---|
| 599 | c
|
|---|
| 600 | c
|
|---|
| 601 | c Geometry:
|
|---|
| 602 | c
|
|---|
| 603 | c
|
|---|
| 604 | c dx0 dx1 dx2 dx3 dx5 dx5 dx6
|
|---|
| 605 | c
|
|---|
| 606 | c bl |<--ll-->|<------lm------>|<---lr--->| br
|
|---|
| 607 | c 0--------x-----|--x---x--------x---x--|-------x-----------0
|
|---|
| 608 | c
|
|---|
| 609 | c left elem. middle elem. right elem.
|
|---|
| 610 | c -1 +1
|
|---|
| 611 | c
|
|---|
| 612 | c
|
|---|
| 613 | c "bl" = (extrapolated) location of Gauss point lx2-1 in left elem.
|
|---|
| 614 | c
|
|---|
| 615 | c "br" = (extrapolated) location of Gauss point 2 in right elem.
|
|---|
| 616 | c
|
|---|
| 617 | c Overlapping Schwarz applied with homogeneous Dirichlet boundary
|
|---|
| 618 | c conditions at "bl" and "br", and with a single d.o.f. extending
|
|---|
| 619 | c in to each adjacent domain.
|
|---|
| 620 | c
|
|---|
| 621 | eps = 1.e-5
|
|---|
| 622 | call rone(dx,nz+3)
|
|---|
| 623 | c
|
|---|
| 624 | c Middle
|
|---|
| 625 | scale = lm/(z(nz)-z(1))
|
|---|
| 626 | do i=1,nz-1
|
|---|
| 627 | dx(i+1) = (z(i+1)-z(i))*scale
|
|---|
| 628 | enddo
|
|---|
| 629 |
|
|---|
| 630 | c Left end
|
|---|
| 631 | if (lbc.eq.0) then
|
|---|
| 632 | dzm0 = z(1) + 1.
|
|---|
| 633 | dxm0 = scale*dzm0
|
|---|
| 634 | dxln = ll - dxm0
|
|---|
| 635 | scalel = dxln/dzm0
|
|---|
| 636 | dx(0) = scalel*(z(2)-z(1))
|
|---|
| 637 | dx(1) = ll
|
|---|
| 638 | elseif (lbc.eq.1) then
|
|---|
| 639 | dx(1) = ll
|
|---|
| 640 | endif
|
|---|
| 641 | c
|
|---|
| 642 | c Right end
|
|---|
| 643 | if (rbc.eq.0) then
|
|---|
| 644 | dzm0 = z(1) + 1.
|
|---|
| 645 | dxm0 = scale*dzm0
|
|---|
| 646 | dxr0 = lr - dxm0
|
|---|
| 647 | scaler = dxr0/dzm0
|
|---|
| 648 | dx(nz+1) = lr
|
|---|
| 649 | dx(nz+2) = scaler*(z(2)-z(1))
|
|---|
| 650 | elseif (rbc.eq.1) then
|
|---|
| 651 | dx(nz+1) = lr
|
|---|
| 652 | endif
|
|---|
| 653 | c
|
|---|
| 654 | return
|
|---|
| 655 | end
|
|---|
| 656 | c-----------------------------------------------------------------------
|
|---|
| 657 | subroutine gen_eigs_A_fem(sf,sft,atd,n,l,lbc,rbc)
|
|---|
| 658 | c
|
|---|
| 659 | c Set up Fast diagonalization solver for FEM on mesh 2
|
|---|
| 660 | c
|
|---|
| 661 | real sf(n,n),sft(n,n),atd(1),l(0:1)
|
|---|
| 662 | integer lbc,rbc
|
|---|
| 663 | c
|
|---|
| 664 | parameter (m=100)
|
|---|
| 665 | real atu(m),ad(m),au(m),c(m),bh(m),li(0:m)
|
|---|
| 666 | c
|
|---|
| 667 | if (n.gt.m) then
|
|---|
| 668 | write(6,*) 'ABORT. Error in gen_eigs_A_fem. Increase m to',n
|
|---|
| 669 | call exitt
|
|---|
| 670 | endif
|
|---|
| 671 | c
|
|---|
| 672 | c Get delta x's
|
|---|
| 673 | c
|
|---|
| 674 | do i=0,n
|
|---|
| 675 | li(i) = 1.0/l(i)
|
|---|
| 676 | enddo
|
|---|
| 677 | c ^ ^
|
|---|
| 678 | c Fill initial arrays, A & B:
|
|---|
| 679 | c
|
|---|
| 680 | call rzero(ad,n)
|
|---|
| 681 | call rzero(au,n-1)
|
|---|
| 682 | call rzero(bh,n)
|
|---|
| 683 | c
|
|---|
| 684 | ie1 = lbc
|
|---|
| 685 | ien = n-rbc
|
|---|
| 686 | do ie=ie1,ien
|
|---|
| 687 | c
|
|---|
| 688 | c il,ir are the left and right endpts of element ie.
|
|---|
| 689 | il = ie
|
|---|
| 690 | ir = ie+1
|
|---|
| 691 | c
|
|---|
| 692 | if (ie.gt.0) ad(il) = ad(il) + li(ie)
|
|---|
| 693 | if (ie.lt.n) ad(ir) = ad(ir) + li(ie)
|
|---|
| 694 | if (ie.gt.0) au(il) = - li(ie)
|
|---|
| 695 | if (ie.gt.0) bh(il) = bh(il) + 0.5 * l(ie)
|
|---|
| 696 | if (ie.lt.n) bh(ir) = bh(ir) + 0.5 * l(ie)
|
|---|
| 697 | enddo
|
|---|
| 698 | c
|
|---|
| 699 | c Take care of bc's (using blasting)
|
|---|
| 700 | bhm = vlmax(bh(2),n-2)/(n-2)
|
|---|
| 701 | ahm = vlmax(ad(2),n-2)/(n-2)
|
|---|
| 702 | c
|
|---|
| 703 | if (lbc.gt.0) then
|
|---|
| 704 | au(1) = 0.
|
|---|
| 705 | ad(1) = ahm
|
|---|
| 706 | bh(1) = bhm
|
|---|
| 707 | endif
|
|---|
| 708 | c
|
|---|
| 709 | if (rbc.gt.0) then
|
|---|
| 710 | au(n-1) = 0.
|
|---|
| 711 | ad(n ) = ahm
|
|---|
| 712 | bh(n ) = bhm
|
|---|
| 713 | endif
|
|---|
| 714 | c
|
|---|
| 715 | c
|
|---|
| 716 | do i=1,n
|
|---|
| 717 | c(i) = sqrt(1.0/bh(i))
|
|---|
| 718 | enddo
|
|---|
| 719 | c ~
|
|---|
| 720 | c Scale rows and columns of A by C: A = CAC
|
|---|
| 721 | c
|
|---|
| 722 | do i=1,n
|
|---|
| 723 | atd(i) = c(i)*ad(i)*c(i)
|
|---|
| 724 | enddo
|
|---|
| 725 | c
|
|---|
| 726 | c Scale upper diagonal
|
|---|
| 727 | c
|
|---|
| 728 | atu(1) = 0.
|
|---|
| 729 | do i=1,n-1
|
|---|
| 730 | atu(i) = c(i)*au(i)*c(i+1)
|
|---|
| 731 | enddo
|
|---|
| 732 | c ~
|
|---|
| 733 | c Compute eigenvalues and eigenvectors of A
|
|---|
| 734 | c
|
|---|
| 735 | call calcz(atd,atu,n,dmax,dmin,sf,ierr)
|
|---|
| 736 | if (ierr.eq.1) then
|
|---|
| 737 | nid = mynode()
|
|---|
| 738 | write(6,6) nid,' czfail:',(l(k),k=0,n)
|
|---|
| 739 | 6 format(i5,a8,1p16e10.2)
|
|---|
| 740 | call exitt
|
|---|
| 741 | endif
|
|---|
| 742 | c
|
|---|
| 743 | c Sort eigenvalues and vectors
|
|---|
| 744 | c
|
|---|
| 745 | call sort(atd,atu,n)
|
|---|
| 746 | call transpose(sft,n,sf,n)
|
|---|
| 747 | do j=1,n
|
|---|
| 748 | call swap(sft(1,j),atu,n,au)
|
|---|
| 749 | enddo
|
|---|
| 750 | call transpose(sf,n,sft,n)
|
|---|
| 751 | c
|
|---|
| 752 | c Make "like" eigenvectors of same sign (for ease of diagnostics)
|
|---|
| 753 | c
|
|---|
| 754 | do j=1,n
|
|---|
| 755 | avg = vlsum(sf(1,j),n)
|
|---|
| 756 | if (avg.lt.0) call chsign(sf(1,j),n)
|
|---|
| 757 | enddo
|
|---|
| 758 | c
|
|---|
| 759 | c Clean up zero eigenvalue
|
|---|
| 760 | c
|
|---|
| 761 | eps = 1.0e-6*dmax
|
|---|
| 762 | do i=1,n
|
|---|
| 763 | c if (atd(i).lt.eps)
|
|---|
| 764 | c $ write(6,*) 'Reset Eig in gen_Afem:',i,n,atd(i)
|
|---|
| 765 | if (atd(i).lt.eps) atd(i) = 0.0
|
|---|
| 766 | enddo
|
|---|
| 767 | c
|
|---|
| 768 | c scale eigenvectors by C:
|
|---|
| 769 | c
|
|---|
| 770 | do i=1,n
|
|---|
| 771 | do j=1,n
|
|---|
| 772 | sf(i,j) = sf(i,j)*c(i)
|
|---|
| 773 | enddo
|
|---|
| 774 | enddo
|
|---|
| 775 | c ^
|
|---|
| 776 | c Orthonormalize eigenvectors w.r.t. B inner-product
|
|---|
| 777 | c
|
|---|
| 778 | do j=1,n
|
|---|
| 779 | alpha = vlsc3(bh,sf(1,j),sf(1,j),n)
|
|---|
| 780 | alpha = 1.0/sqrt(alpha)
|
|---|
| 781 | call cmult(sf(1,j),alpha,n)
|
|---|
| 782 | enddo
|
|---|
| 783 | c
|
|---|
| 784 | c Diagnostics
|
|---|
| 785 | c
|
|---|
| 786 | c do j=1,n
|
|---|
| 787 | c do i=1,n
|
|---|
| 788 | c sft(i,j) = vlsc3(bh,sf(1,i),sf(1,j),n)
|
|---|
| 789 | c enddo
|
|---|
| 790 | c enddo
|
|---|
| 791 | c
|
|---|
| 792 | c n8 = min(n,8)
|
|---|
| 793 | c do i=1,n
|
|---|
| 794 | c write(6,2) (sft(i,j),j=1,n8)
|
|---|
| 795 | c enddo
|
|---|
| 796 | c 2 format('Id:',1p8e12.4)
|
|---|
| 797 | c
|
|---|
| 798 | call transpose(sft,n,sf,n)
|
|---|
| 799 | return
|
|---|
| 800 | end
|
|---|
| 801 | c-----------------------------------------------------------------------
|
|---|
| 802 | subroutine get_fast_bc(lbr,rbr,lbs,rbs,lbt,rbt,e,bsym,ierr)
|
|---|
| 803 | integer lbr,rbr,lbs,rbs,lbt,rbt,e,bsym
|
|---|
| 804 | c
|
|---|
| 805 | include 'SIZE'
|
|---|
| 806 | include 'INPUT'
|
|---|
| 807 | include 'PARALLEL'
|
|---|
| 808 | include 'TOPOL'
|
|---|
| 809 | include 'TSTEP'
|
|---|
| 810 | c
|
|---|
| 811 | integer fbc(6)
|
|---|
| 812 | c
|
|---|
| 813 | c ibc = 0 <==> Dirichlet
|
|---|
| 814 | c ibc = 1 <==> Dirichlet, outflow (no extension)
|
|---|
| 815 | c ibc = 2 <==> Neumann,
|
|---|
| 816 |
|
|---|
| 817 |
|
|---|
| 818 | do iface=1,2*ldim
|
|---|
| 819 | ied = eface(iface)
|
|---|
| 820 | ibc = -1
|
|---|
| 821 |
|
|---|
| 822 | if (ifmhd) call mhd_bc_dn(ibc,iface,e) ! can be overwritten by 'mvn'
|
|---|
| 823 |
|
|---|
| 824 | if (cbc(ied,e,ifield).eq.' ') ibc = 0
|
|---|
| 825 | if (cbc(ied,e,ifield).eq.'E ') ibc = 0
|
|---|
| 826 | if (cbc(ied,e,ifield).eq.'msi') ibc = 0
|
|---|
| 827 | if (cbc(ied,e,ifield).eq.'MSI') ibc = 0
|
|---|
| 828 | if (cbc(ied,e,ifield).eq.'P ') ibc = 0
|
|---|
| 829 | if (cbc(ied,e,ifield).eq.'p ') ibc = 0
|
|---|
| 830 | if (cbc(ied,e,ifield).eq.'O ') ibc = 1
|
|---|
| 831 | if (cbc(ied,e,ifield).eq.'ON ') ibc = 1
|
|---|
| 832 | if (cbc(ied,e,ifield).eq.'o ') ibc = 1
|
|---|
| 833 | if (cbc(ied,e,ifield).eq.'on ') ibc = 1
|
|---|
| 834 | if (cbc(ied,e,ifield).eq.'MS ') ibc = 1
|
|---|
| 835 | if (cbc(ied,e,ifield).eq.'ms ') ibc = 1
|
|---|
| 836 | if (cbc(ied,e,ifield).eq.'MM ') ibc = 1
|
|---|
| 837 | if (cbc(ied,e,ifield).eq.'mm ') ibc = 1
|
|---|
| 838 | if (cbc(ied,e,ifield).eq.'mv ') ibc = 2
|
|---|
| 839 | if (cbc(ied,e,ifield).eq.'mvn') ibc = 2
|
|---|
| 840 | if (cbc(ied,e,ifield).eq.'v ') ibc = 2
|
|---|
| 841 | if (cbc(ied,e,ifield).eq.'V ') ibc = 2
|
|---|
| 842 | if (cbc(ied,e,ifield).eq.'W ') ibc = 2
|
|---|
| 843 | if (cbc(ied,e,ifield).eq.'SYM') ibc = bsym
|
|---|
| 844 | if (cbc(ied,e,ifield).eq.'SL ') ibc = 2
|
|---|
| 845 | if (cbc(ied,e,ifield).eq.'sl ') ibc = 2
|
|---|
| 846 | if (cbc(ied,e,ifield).eq.'SHL') ibc = 2
|
|---|
| 847 | if (cbc(ied,e,ifield).eq.'shl') ibc = 2
|
|---|
| 848 | if (cbc(ied,e,ifield).eq.'A ') ibc = 2
|
|---|
| 849 | if (cbc(ied,e,ifield).eq.'S ') ibc = 2
|
|---|
| 850 | if (cbc(ied,e,ifield).eq.'s ') ibc = 2
|
|---|
| 851 | if (cbc(ied,e,ifield).eq.'J ') ibc = 0
|
|---|
| 852 | if (cbc(ied,e,ifield).eq.'SP ') ibc = 0
|
|---|
| 853 |
|
|---|
| 854 | fbc(iface) = ibc
|
|---|
| 855 |
|
|---|
| 856 | if (ierr.eq.-1) write(6,1) ibc,ied,e,ifield,cbc(ied,e,ifield)
|
|---|
| 857 | 1 format(2i3,i8,i3,2x,a3,' get_fast_bc_error')
|
|---|
| 858 |
|
|---|
| 859 | enddo
|
|---|
| 860 |
|
|---|
| 861 | if (ierr.eq.-1) call exitti('Error A get_fast_bc$',e)
|
|---|
| 862 |
|
|---|
| 863 | lbr = fbc(1)
|
|---|
| 864 | rbr = fbc(2)
|
|---|
| 865 | lbs = fbc(3)
|
|---|
| 866 | rbs = fbc(4)
|
|---|
| 867 | lbt = fbc(5)
|
|---|
| 868 | rbt = fbc(6)
|
|---|
| 869 |
|
|---|
| 870 | ierr = 0
|
|---|
| 871 | if (ibc.lt.0) ierr = lglel(e)
|
|---|
| 872 |
|
|---|
| 873 | c write(6,6) e,lbr,rbr,lbs,rbs,(cbc(k,e,ifield),k=1,4)
|
|---|
| 874 | c 6 format(i5,2x,4i3,3x,4(1x,a3),' get_fast_bc')
|
|---|
| 875 |
|
|---|
| 876 | return
|
|---|
| 877 | end
|
|---|
| 878 | c-----------------------------------------------------------------------
|
|---|
| 879 | subroutine outv(x,n,name3)
|
|---|
| 880 | character*3 name3
|
|---|
| 881 | real x(1)
|
|---|
| 882 | c
|
|---|
| 883 | nn = min (n,10)
|
|---|
| 884 | write(6,6) name3,(x(i),i=1,nn)
|
|---|
| 885 | 6 format(a3,10f12.6)
|
|---|
| 886 | c
|
|---|
| 887 | return
|
|---|
| 888 | end
|
|---|
| 889 | c-----------------------------------------------------------------------
|
|---|
| 890 | subroutine outmat(a,m,n,name6,ie)
|
|---|
| 891 | real a(m,n)
|
|---|
| 892 | character*6 name6
|
|---|
| 893 | c
|
|---|
| 894 | write(6,*)
|
|---|
| 895 | write(6,*) ie,' matrix: ',name6,m,n
|
|---|
| 896 | n12 = min(n,12)
|
|---|
| 897 | do i=1,m
|
|---|
| 898 | write(6,6) ie,name6,(a(i,j),j=1,n12)
|
|---|
| 899 | enddo
|
|---|
| 900 | 6 format(i3,1x,a6,12f9.5)
|
|---|
| 901 | write(6,*)
|
|---|
| 902 | return
|
|---|
| 903 | end
|
|---|
| 904 | c-----------------------------------------------------------------------
|
|---|
| 905 | subroutine set_up_fast_1D_fem_ax
|
|---|
| 906 | $ (s,lam,n,lbc,rbc,ll,lm,lr,z,y,nz,ie)
|
|---|
| 907 | real s(1),lam(1),ll,lm,lr,z(1),y(1)
|
|---|
| 908 | integer lbc,rbc
|
|---|
| 909 | c
|
|---|
| 910 | parameter (m=100)
|
|---|
| 911 | real dx(0:m)
|
|---|
| 912 | integer icalld
|
|---|
| 913 | save icalld
|
|---|
| 914 | data icalld/0/
|
|---|
| 915 | c
|
|---|
| 916 | icalld=icalld+1
|
|---|
| 917 | c
|
|---|
| 918 | if (nz.gt.m-3) then
|
|---|
| 919 | write(6,*) 'ABORT. Error in set_up_fast_1D_fem. Increase m to'
|
|---|
| 920 | $ , nz
|
|---|
| 921 | call exitt
|
|---|
| 922 | endif
|
|---|
| 923 | c
|
|---|
| 924 | c In the present scheme, each element is viewed as a d-fold
|
|---|
| 925 | c tensor of (1+nz+1) arrays, even if funky bc's are applied
|
|---|
| 926 | c on either end of the 1D array.
|
|---|
| 927 | c
|
|---|
| 928 | n = nz+2
|
|---|
| 929 | c
|
|---|
| 930 | c Compute spacing, dx()
|
|---|
| 931 | c
|
|---|
| 932 | call set_up_1D_geom_ax(dx,lbc,rbc,ll,lm,lr,z,y,nz)
|
|---|
| 933 | c
|
|---|
| 934 | nn1 = n*n + 1
|
|---|
| 935 | call gen_eigs_A_fem_ax(s,s(nn1),lam,n,dx,y,lbc,rbc)
|
|---|
| 936 | c
|
|---|
| 937 | return
|
|---|
| 938 | end
|
|---|
| 939 | c-----------------------------------------------------------------------
|
|---|
| 940 | subroutine set_up_1D_geom_ax(dx,lbc,rbc,ll,lm,lr,z,y,nz)
|
|---|
| 941 | c
|
|---|
| 942 | real dx(0:1),ll,lm,lr,z(2),y(1)
|
|---|
| 943 | integer lbc,rbc
|
|---|
| 944 | c
|
|---|
| 945 | c
|
|---|
| 946 | c Set up the 1D geometry for the tensor-product based overlapping Schwarz
|
|---|
| 947 | c
|
|---|
| 948 | c Upon return:
|
|---|
| 949 | c
|
|---|
| 950 | c dx() contains the spacing required to set up the stiffness matrix.
|
|---|
| 951 | c
|
|---|
| 952 | c
|
|---|
| 953 | c Input:
|
|---|
| 954 | c
|
|---|
| 955 | c lbc (rbc) is 0 if the left (right) BC is Dirichlet, 1 if Neumann.
|
|---|
| 956 | c
|
|---|
| 957 | c ll is the space between the left-most Gauss point of the middle
|
|---|
| 958 | c element and the right-most Gauss point of the LEFT element
|
|---|
| 959 | c
|
|---|
| 960 | c lm is the space between the left-most Gauss point of the middle
|
|---|
| 961 | c element and the right-most Gauss point of the MIDDLE element
|
|---|
| 962 | c
|
|---|
| 963 | c lr is the space between the right-most Gauss point of the middle
|
|---|
| 964 | c element and the left-most Gauss point of the RIGHT element
|
|---|
| 965 | c
|
|---|
| 966 | c --- if ll (lr) is very small (0), it indicates that there is no
|
|---|
| 967 | c left (right) spacing, and that the left (right) BC is Neumann.
|
|---|
| 968 | c
|
|---|
| 969 | c
|
|---|
| 970 | c z() is the array of nz Gauss points on the interval ]-1,1[.
|
|---|
| 971 | c
|
|---|
| 972 | c Boundary conditions:
|
|---|
| 973 | c
|
|---|
| 974 | c bc = 0 -- std. Dirichlet bc applied 2 points away from interior
|
|---|
| 975 | c bc = 1 -- Dirichlet bc applied 1 point away from interior (outflow)
|
|---|
| 976 | c bc = 2 -- Neumann bc applied on interior point (W,v,V,SYM,...)
|
|---|
| 977 | c
|
|---|
| 978 | c
|
|---|
| 979 | c
|
|---|
| 980 | c Geometry:
|
|---|
| 981 | c
|
|---|
| 982 | c
|
|---|
| 983 | c dx0 dx1 dx2 dx3 dx5 dx5 dx6
|
|---|
| 984 | c
|
|---|
| 985 | c bl |<--ll-->|<------lm------>|<---lr--->| br
|
|---|
| 986 | c 0--------x-----|--x---x--------x---x--|-------x-----------0
|
|---|
| 987 | c
|
|---|
| 988 | c left elem. middle elem. right elem.
|
|---|
| 989 | c -1 +1
|
|---|
| 990 | c
|
|---|
| 991 | c
|
|---|
| 992 | c "bl" = (extrapolated) location of Gauss point lx2-1 in left elem.
|
|---|
| 993 | c
|
|---|
| 994 | c "br" = (extrapolated) location of Gauss point 2 in right elem.
|
|---|
| 995 | c
|
|---|
| 996 | c Overlapping Schwarz applied with homogeneous Dirichlet boundary
|
|---|
| 997 | c conditions at "bl" and "br", and with a single d.o.f. extending
|
|---|
| 998 | c in to each adjacent domain.
|
|---|
| 999 | c
|
|---|
| 1000 | eps = 1.e-5
|
|---|
| 1001 | call rone(dx,nz+3)
|
|---|
| 1002 | c
|
|---|
| 1003 | c Middle
|
|---|
| 1004 | scale = lm/(z(nz)-z(1))
|
|---|
| 1005 | do i=1,nz-1
|
|---|
| 1006 | dx(i+1) = (z(i+1)-z(i))*scale
|
|---|
| 1007 | enddo
|
|---|
| 1008 |
|
|---|
| 1009 | c Left end
|
|---|
| 1010 | if (lbc.eq.0) then
|
|---|
| 1011 | dzm0 = z(1) + 1.
|
|---|
| 1012 | dxm0 = scale*dzm0
|
|---|
| 1013 | dxln = ll - dxm0
|
|---|
| 1014 | scalel = dxln/dzm0
|
|---|
| 1015 | dx(0) = scalel*(z(2)-z(1))
|
|---|
| 1016 | dx(1) = ll
|
|---|
| 1017 | elseif (lbc.eq.1) then
|
|---|
| 1018 | dx(1) = ll
|
|---|
| 1019 | endif
|
|---|
| 1020 | c
|
|---|
| 1021 | c Right end
|
|---|
| 1022 | if (rbc.eq.0) then
|
|---|
| 1023 | dzm0 = z(1) + 1.
|
|---|
| 1024 | dxm0 = scale*dzm0
|
|---|
| 1025 | dxr0 = lr - dxm0
|
|---|
| 1026 | scaler = dxr0/dzm0
|
|---|
| 1027 | dx(nz+1) = lr
|
|---|
| 1028 | dx(nz+2) = scaler*(z(2)-z(1))
|
|---|
| 1029 | elseif (rbc.eq.1) then
|
|---|
| 1030 | dx(nz+1) = lr
|
|---|
| 1031 | endif
|
|---|
| 1032 | c
|
|---|
| 1033 | return
|
|---|
| 1034 | end
|
|---|
| 1035 | c-----------------------------------------------------------------------
|
|---|
| 1036 | subroutine gen_eigs_A_fem_ax(sf,sft,atd,n,l,y,lbc,rbc)
|
|---|
| 1037 | c
|
|---|
| 1038 | c Set up Fast diagonalization solver for FEM on mesh 2
|
|---|
| 1039 | c
|
|---|
| 1040 | real sf(n,n),sft(n,n),atd(1),l(0:1),y(1)
|
|---|
| 1041 | integer lbc,rbc
|
|---|
| 1042 | c
|
|---|
| 1043 | parameter (m=100)
|
|---|
| 1044 | real atu(m),ad(m),au(m),c(m),bh(m),li(0:m)
|
|---|
| 1045 | c
|
|---|
| 1046 | if (n.gt.m) then
|
|---|
| 1047 | write(6,*) 'ABORT. Error in gen_eigs_A_fem. Increase m to',n
|
|---|
| 1048 | call exitt
|
|---|
| 1049 | endif
|
|---|
| 1050 | c
|
|---|
| 1051 | c Get delta x's
|
|---|
| 1052 | c
|
|---|
| 1053 | do i=0,n
|
|---|
| 1054 | li(i) = 1.0/l(i)
|
|---|
| 1055 | enddo
|
|---|
| 1056 | c ^ ^
|
|---|
| 1057 | c Fill initial arrays, A & B:
|
|---|
| 1058 | c
|
|---|
| 1059 | call rzero(ad,n)
|
|---|
| 1060 | call rzero(au,n-1)
|
|---|
| 1061 | call rzero(bh,n)
|
|---|
| 1062 | c
|
|---|
| 1063 | ie1 = lbc
|
|---|
| 1064 | ien = n-rbc
|
|---|
| 1065 | do ie=ie1,ien
|
|---|
| 1066 | c
|
|---|
| 1067 | c il,ir are the left and right endpts of element ie.
|
|---|
| 1068 | il = ie
|
|---|
| 1069 | ir = ie+1
|
|---|
| 1070 | c
|
|---|
| 1071 | if (ie.gt.0) ad(il) = ad(il) + li(ie)
|
|---|
| 1072 | if (ie.lt.n) ad(ir) = ad(ir) + li(ie)
|
|---|
| 1073 | if (ie.gt.0) au(il) = - li(ie)
|
|---|
| 1074 | if (ie.gt.0) bh(il) = bh(il) + 0.5 * l(ie)
|
|---|
| 1075 | if (ie.lt.n) bh(ir) = bh(ir) + 0.5 * l(ie)
|
|---|
| 1076 | enddo
|
|---|
| 1077 | c
|
|---|
| 1078 | c Take care of bc's (using blasting)
|
|---|
| 1079 | bhm = vlmax(bh(2),n-2)/(n-2)
|
|---|
| 1080 | ahm = vlmax(ad(2),n-2)/(n-2)
|
|---|
| 1081 | c
|
|---|
| 1082 | if (lbc.gt.0) then
|
|---|
| 1083 | au(1) = 0.
|
|---|
| 1084 | ad(1) = ahm
|
|---|
| 1085 | bh(1) = bhm
|
|---|
| 1086 | endif
|
|---|
| 1087 | c
|
|---|
| 1088 | if (rbc.gt.0) then
|
|---|
| 1089 | au(n-1) = 0.
|
|---|
| 1090 | ad(n ) = ahm
|
|---|
| 1091 | bh(n ) = bhm
|
|---|
| 1092 | endif
|
|---|
| 1093 | c
|
|---|
| 1094 | c
|
|---|
| 1095 | do i=1,n
|
|---|
| 1096 | c(i) = sqrt(1.0/bh(i))
|
|---|
| 1097 | enddo
|
|---|
| 1098 | c ~
|
|---|
| 1099 | c Scale rows and columns of A by C: A = CAC
|
|---|
| 1100 | c
|
|---|
| 1101 | do i=1,n
|
|---|
| 1102 | atd(i) = c(i)*ad(i)*c(i)
|
|---|
| 1103 | enddo
|
|---|
| 1104 | c
|
|---|
| 1105 | c Scale upper diagonal
|
|---|
| 1106 | c
|
|---|
| 1107 | atu(1) = 0.
|
|---|
| 1108 | do i=1,n-1
|
|---|
| 1109 | atu(i) = c(i)*au(i)*c(i+1)
|
|---|
| 1110 | enddo
|
|---|
| 1111 | c ~
|
|---|
| 1112 | c Compute eigenvalues and eigenvectors of A
|
|---|
| 1113 | c
|
|---|
| 1114 | call calcz(atd,atu,n,dmax,dmin,sf,ierr)
|
|---|
| 1115 | if (ierr.eq.1) then
|
|---|
| 1116 | nid = mynode()
|
|---|
| 1117 | write(6,6) nid,' czfail2:',(l(k),k=0,n)
|
|---|
| 1118 | 6 format(i5,a8,1p16e10.2)
|
|---|
| 1119 | call exitt
|
|---|
| 1120 | endif
|
|---|
| 1121 | c
|
|---|
| 1122 | c Sort eigenvalues and vectors
|
|---|
| 1123 | c
|
|---|
| 1124 | call sort(atd,atu,n)
|
|---|
| 1125 | call transpose(sft,n,sf,n)
|
|---|
| 1126 | do j=1,n
|
|---|
| 1127 | call swap(sft(1,j),atu,n,au)
|
|---|
| 1128 | enddo
|
|---|
| 1129 | call transpose(sf,n,sft,n)
|
|---|
| 1130 | c
|
|---|
| 1131 | c Make "like" eigenvectors of same sign (for ease of diagnostics)
|
|---|
| 1132 | c
|
|---|
| 1133 | do j=1,n
|
|---|
| 1134 | avg = vlsum(sf(1,j),n)
|
|---|
| 1135 | if (avg.lt.0) call chsign(sf(1,j),n)
|
|---|
| 1136 | enddo
|
|---|
| 1137 | c
|
|---|
| 1138 | c Clean up zero eigenvalue
|
|---|
| 1139 | c
|
|---|
| 1140 | eps = 1.0e-6*dmax
|
|---|
| 1141 | do i=1,n
|
|---|
| 1142 | c if (atd(i).lt.eps)
|
|---|
| 1143 | c $ write(6,*) 'Reset Eig in gen_Afm_ax:',i,n,atd(i)
|
|---|
| 1144 | if (atd(i).lt.eps) atd(i) = 0.0
|
|---|
| 1145 | enddo
|
|---|
| 1146 | c
|
|---|
| 1147 | c scale eigenvectors by C:
|
|---|
| 1148 | c
|
|---|
| 1149 | do i=1,n
|
|---|
| 1150 | do j=1,n
|
|---|
| 1151 | sf(i,j) = sf(i,j)*c(i)
|
|---|
| 1152 | enddo
|
|---|
| 1153 | enddo
|
|---|
| 1154 | c ^
|
|---|
| 1155 | c Orthonormalize eigenvectors w.r.t. B inner-product
|
|---|
| 1156 | c
|
|---|
| 1157 | do j=1,n
|
|---|
| 1158 | alpha = vlsc3(bh,sf(1,j),sf(1,j),n)
|
|---|
| 1159 | alpha = 1.0/sqrt(alpha)
|
|---|
| 1160 | call cmult(sf(1,j),alpha,n)
|
|---|
| 1161 | enddo
|
|---|
| 1162 | c
|
|---|
| 1163 | c Diagnostics
|
|---|
| 1164 | c
|
|---|
| 1165 | c do j=1,n
|
|---|
| 1166 | c do i=1,n
|
|---|
| 1167 | c sft(i,j) = vlsc3(bh,sf(1,i),sf(1,j),n)
|
|---|
| 1168 | c enddo
|
|---|
| 1169 | c enddo
|
|---|
| 1170 | c
|
|---|
| 1171 | c n8 = min(n,8)
|
|---|
| 1172 | c do i=1,n
|
|---|
| 1173 | c write(6,2) (sft(i,j),j=1,n8)
|
|---|
| 1174 | c enddo
|
|---|
| 1175 | c 2 format('Id:',1p8e12.4)
|
|---|
| 1176 | c
|
|---|
| 1177 | call transpose(sft,n,sf,n)
|
|---|
| 1178 | return
|
|---|
| 1179 | end
|
|---|
| 1180 | c-----------------------------------------------------------------------
|
|---|
| 1181 | subroutine load_semhat_weighted ! Fills the SEMHAT arrays
|
|---|
| 1182 | c
|
|---|
| 1183 | c Note that this routine performs the following matrix multiplies
|
|---|
| 1184 | c after getting the matrices back from semhat:
|
|---|
| 1185 | c
|
|---|
| 1186 | c dgl = bgl dgl
|
|---|
| 1187 | c jgl = bgl jgl
|
|---|
| 1188 | c
|
|---|
| 1189 | include 'SIZE'
|
|---|
| 1190 | include 'SEMHAT'
|
|---|
| 1191 | c
|
|---|
| 1192 | nr = lx1-1
|
|---|
| 1193 | call semhat(ah,bh,ch,dh,zh,dph,jph,bgl,zglhat,dgl,jgl,nr,wh)
|
|---|
| 1194 | call do_semhat_weight(jgl,dgl,bgl,nr)
|
|---|
| 1195 | c
|
|---|
| 1196 | return
|
|---|
| 1197 | end
|
|---|
| 1198 | c----------------------------------------------------------------------
|
|---|
| 1199 | subroutine do_semhat_weight(jgl,dgl,bgl,n)
|
|---|
| 1200 | real bgl(1:n-1),jgl(1:n-1,0:n),dgl(1:n-1,0:n)
|
|---|
| 1201 |
|
|---|
| 1202 | do j=0,n
|
|---|
| 1203 | do i=1,n-1
|
|---|
| 1204 | jgl(i,j)=bgl(i)*jgl(i,j)
|
|---|
| 1205 | enddo
|
|---|
| 1206 | enddo
|
|---|
| 1207 | do j=0,n
|
|---|
| 1208 | do i=1,n-1
|
|---|
| 1209 | dgl(i,j)=bgl(i)*dgl(i,j)
|
|---|
| 1210 | enddo
|
|---|
| 1211 | enddo
|
|---|
| 1212 | return
|
|---|
| 1213 | end
|
|---|
| 1214 | c-----------------------------------------------------------------------
|
|---|
| 1215 | subroutine semhat(a,b,c,d,z,dgll,jgll,bgl,zgl,dgl,jgl,n,w)
|
|---|
| 1216 | c
|
|---|
| 1217 | c Generate matrices for single element, 1D operators:
|
|---|
| 1218 | c
|
|---|
| 1219 | c a = Laplacian
|
|---|
| 1220 | c b = diagonal mass matrix
|
|---|
| 1221 | c c = convection operator b*d
|
|---|
| 1222 | c d = derivative matrix
|
|---|
| 1223 | c dgll = derivative matrix, mapping from pressure nodes to velocity
|
|---|
| 1224 | c jgll = interpolation matrix, mapping from pressure nodes to velocity
|
|---|
| 1225 | c z = GLL points
|
|---|
| 1226 | c
|
|---|
| 1227 | c zgl = GL points
|
|---|
| 1228 | c bgl = diagonal mass matrix on GL
|
|---|
| 1229 | c dgl = derivative matrix, mapping from velocity nodes to pressure
|
|---|
| 1230 | c jgl = interpolation matrix, mapping from velocity nodes to pressure
|
|---|
| 1231 | c
|
|---|
| 1232 | c n = polynomial degree (velocity space)
|
|---|
| 1233 | c w = work array of size 2*n+2
|
|---|
| 1234 | c
|
|---|
| 1235 | c Currently, this is set up for pressure nodes on the interior GLL pts.
|
|---|
| 1236 | c
|
|---|
| 1237 | c
|
|---|
| 1238 | real a(0:n,0:n),b(0:n),c(0:n,0:n),d(0:n,0:n),z(0:n)
|
|---|
| 1239 | real dgll(0:n,1:n-1),jgll(0:n,1:n-1)
|
|---|
| 1240 | c
|
|---|
| 1241 | real bgl(1:n-1),zgl(1:n-1)
|
|---|
| 1242 | real dgl(1:n-1,0:n),jgl(1:n-1,0:n)
|
|---|
| 1243 | c
|
|---|
| 1244 | real w(0:2*n+1)
|
|---|
| 1245 | c
|
|---|
| 1246 | np = n+1
|
|---|
| 1247 | nm = n-1
|
|---|
| 1248 | n2 = n-2
|
|---|
| 1249 | c
|
|---|
| 1250 | call zwgll (z,b,np)
|
|---|
| 1251 | c
|
|---|
| 1252 | do i=0,n
|
|---|
| 1253 | call fd_weights_full(z(i),z,n,1,w)
|
|---|
| 1254 | do j=0,n
|
|---|
| 1255 | d(i,j) = w(j+np) ! Derivative matrix
|
|---|
| 1256 | enddo
|
|---|
| 1257 | enddo
|
|---|
| 1258 |
|
|---|
| 1259 | if (n.eq.1) return ! No interpolation for n=1
|
|---|
| 1260 |
|
|---|
| 1261 | do i=0,n
|
|---|
| 1262 | call fd_weights_full(z(i),z(1),n2,1,w(1))
|
|---|
| 1263 | do j=1,nm
|
|---|
| 1264 | jgll(i,j) = w(j ) ! Interpolation matrix
|
|---|
| 1265 | dgll(i,j) = w(j+nm) ! Derivative matrix
|
|---|
| 1266 | enddo
|
|---|
| 1267 | enddo
|
|---|
| 1268 | c
|
|---|
| 1269 | call rzero(a,np*np)
|
|---|
| 1270 | do j=0,n
|
|---|
| 1271 | do i=0,n
|
|---|
| 1272 | do k=0,n
|
|---|
| 1273 | a(i,j) = a(i,j) + d(k,i)*b(k)*d(k,j)
|
|---|
| 1274 | enddo
|
|---|
| 1275 | c(i,j) = b(i)*d(i,j)
|
|---|
| 1276 | enddo
|
|---|
| 1277 | enddo
|
|---|
| 1278 | c
|
|---|
| 1279 | call zwgl (zgl,bgl,nm)
|
|---|
| 1280 | c
|
|---|
| 1281 | do i=1,n-1
|
|---|
| 1282 | call fd_weights_full(zgl(i),z,n,1,w)
|
|---|
| 1283 | do j=0,n
|
|---|
| 1284 | jgl(i,j) = w(j ) ! Interpolation matrix
|
|---|
| 1285 | dgl(i,j) = w(j+np) ! Derivative matrix
|
|---|
| 1286 | enddo
|
|---|
| 1287 | enddo
|
|---|
| 1288 | c
|
|---|
| 1289 | return
|
|---|
| 1290 | end
|
|---|
| 1291 | c-----------------------------------------------------------------------
|
|---|
| 1292 | subroutine fd_weights_full(xx,x,n,m,c)
|
|---|
| 1293 | c
|
|---|
| 1294 | c This routine evaluates the derivative based on all points
|
|---|
| 1295 | c in the stencils. It is more memory efficient than "fd_weights"
|
|---|
| 1296 | c
|
|---|
| 1297 | c This set of routines comes from the appendix of
|
|---|
| 1298 | c A Practical Guide to Pseudospectral Methods, B. Fornberg
|
|---|
| 1299 | c Cambridge Univ. Press, 1996. (pff)
|
|---|
| 1300 | c
|
|---|
| 1301 | c Input parameters:
|
|---|
| 1302 | c xx -- point at wich the approximations are to be accurate
|
|---|
| 1303 | c x -- array of x-ordinates: x(0:n)
|
|---|
| 1304 | c n -- polynomial degree of interpolant (# of points := n+1)
|
|---|
| 1305 | c m -- highest order of derivative to be approxxmated at xi
|
|---|
| 1306 | c
|
|---|
| 1307 | c Output:
|
|---|
| 1308 | c c -- set of coefficients c(0:n,0:m).
|
|---|
| 1309 | c c(j,k) is to be applied at x(j) when
|
|---|
| 1310 | c the kth derivative is approxxmated by a
|
|---|
| 1311 | c stencil extending over x(0),x(1),...x(n).
|
|---|
| 1312 | c
|
|---|
| 1313 | c
|
|---|
| 1314 | real x(0:n),c(0:n,0:m)
|
|---|
| 1315 | c
|
|---|
| 1316 | c1 = 1.
|
|---|
| 1317 | c4 = x(0) - xx
|
|---|
| 1318 | c
|
|---|
| 1319 | do k=0,m
|
|---|
| 1320 | do j=0,n
|
|---|
| 1321 | c(j,k) = 0.
|
|---|
| 1322 | enddo
|
|---|
| 1323 | enddo
|
|---|
| 1324 | c(0,0) = 1.
|
|---|
| 1325 | c
|
|---|
| 1326 | do i=1,n
|
|---|
| 1327 | mn = min(i,m)
|
|---|
| 1328 | c2 = 1.
|
|---|
| 1329 | c5 = c4
|
|---|
| 1330 | c4 = x(i)-xx
|
|---|
| 1331 | do j=0,i-1
|
|---|
| 1332 | c3 = x(i)-x(j)
|
|---|
| 1333 | c2 = c2*c3
|
|---|
| 1334 | do k=mn,1,-1
|
|---|
| 1335 | c(i,k) = c1*(k*c(i-1,k-1)-c5*c(i-1,k))/c2
|
|---|
| 1336 | enddo
|
|---|
| 1337 | c(i,0) = -c1*c5*c(i-1,0)/c2
|
|---|
| 1338 | do k=mn,1,-1
|
|---|
| 1339 | c(j,k) = (c4*c(j,k)-k*c(j,k-1))/c3
|
|---|
| 1340 | enddo
|
|---|
| 1341 | c(j,0) = c4*c(j,0)/c3
|
|---|
| 1342 | enddo
|
|---|
| 1343 | c1 = c2
|
|---|
| 1344 | enddo
|
|---|
| 1345 | c call outmat(c,n+1,m+1,'fdw',n)
|
|---|
| 1346 | return
|
|---|
| 1347 | end
|
|---|
| 1348 | c-----------------------------------------------------------------------
|
|---|
| 1349 | subroutine set_up_fast_1D_sem(s,lam,n,lbc,rbc,ll,lm,lr,ie)
|
|---|
| 1350 | include 'SIZE'
|
|---|
| 1351 | include 'SEMHAT'
|
|---|
| 1352 | c
|
|---|
| 1353 | common /fast1dsem/ g(lr2),w(lr2)
|
|---|
| 1354 | c
|
|---|
| 1355 | real g,w
|
|---|
| 1356 | real s(1),lam(1),ll,lm,lr
|
|---|
| 1357 | integer lbc,rbc
|
|---|
| 1358 |
|
|---|
| 1359 | integer bb0,bb1,eb0,eb1,n,n1
|
|---|
| 1360 | logical l,r
|
|---|
| 1361 |
|
|---|
| 1362 | n=lx1-1
|
|---|
| 1363 | !bcs on E are from normal vel component
|
|---|
| 1364 | if(lbc.eq.2 .or. lbc.eq.3) then !wall,sym - dirichlet velocity
|
|---|
| 1365 | eb0=1
|
|---|
| 1366 | else !outflow,element - neumann velocity
|
|---|
| 1367 | eb0=0
|
|---|
| 1368 | endif
|
|---|
| 1369 | if(rbc.eq.2 .or. rbc.eq.3) then !wall,sym - dirichlet velocity
|
|---|
| 1370 | eb1=n-1
|
|---|
| 1371 | else !outflow,element - neumann velocity
|
|---|
| 1372 | eb1=n
|
|---|
| 1373 | endif
|
|---|
| 1374 | !bcs on B are from tangent vel component
|
|---|
| 1375 | if(lbc.eq.2) then !wall - dirichlet velocity
|
|---|
| 1376 | bb0=1
|
|---|
| 1377 | else !outflow,element,sym - neumann velocity
|
|---|
| 1378 | bb0=0
|
|---|
| 1379 | endif
|
|---|
| 1380 | if(rbc.eq.2) then !wall - dirichlet velocity
|
|---|
| 1381 | bb1=n-1
|
|---|
| 1382 | else !outflow,element,sym - neumann velocity
|
|---|
| 1383 | bb1=n
|
|---|
| 1384 | endif
|
|---|
| 1385 | c
|
|---|
| 1386 | l = (lbc.eq.0)
|
|---|
| 1387 | r = (rbc.eq.0)
|
|---|
| 1388 | c
|
|---|
| 1389 | c calculate E tilde operator
|
|---|
| 1390 | call set_up_fast_1D_sem_op(s,eb0,eb1,l,r,ll,lm,lr,bh,dgl,0)
|
|---|
| 1391 | c call outmat(s,n+1,n+1,' Et ',ie)
|
|---|
| 1392 | c calculate B tilde operator
|
|---|
| 1393 | call set_up_fast_1D_sem_op(g,bb0,bb1,l,r,ll,lm,lr,bh,jgl,1)
|
|---|
| 1394 | c call outmat(g,n+1,n+1,' Bt ',ie)
|
|---|
| 1395 |
|
|---|
| 1396 | n=n+1
|
|---|
| 1397 | call generalev(s,g,lam,n,w)
|
|---|
| 1398 | if(.not.l) call row_zero(s,n,n,1)
|
|---|
| 1399 | if(.not.r) call row_zero(s,n,n,n)
|
|---|
| 1400 | call transpose(s(n*n+1),n,s,n) ! compute the transpose of s
|
|---|
| 1401 |
|
|---|
| 1402 | c call outmat (s,n,n,' S ',ie)
|
|---|
| 1403 | c call outmat (s(n*n+1),n,n,' St ',1)
|
|---|
| 1404 | c call exitt
|
|---|
| 1405 | return
|
|---|
| 1406 | end
|
|---|
| 1407 | c-----------------------------------------------------------------------
|
|---|
| 1408 | subroutine set_up_fast_1D_sem_op(g,b0,b1,l,r,ll,lm,lr,bh,jgl,jscl)
|
|---|
| 1409 | c -1 T
|
|---|
| 1410 | c G = J B J
|
|---|
| 1411 | c
|
|---|
| 1412 | c gives the inexact restriction of this matrix to
|
|---|
| 1413 | c an element plus one node on either side
|
|---|
| 1414 | c
|
|---|
| 1415 | c g - the output matrix
|
|---|
| 1416 | c b0, b1 - the range for Bhat indices for the element
|
|---|
| 1417 | c (enforces boundary conditions)
|
|---|
| 1418 | c l, r - whether there is a left or right neighbor
|
|---|
| 1419 | c ll,lm,lr - lengths of left, middle, and right elements
|
|---|
| 1420 | c bh - hat matrix for B
|
|---|
| 1421 | c jgl - hat matrix for J (should map vel to pressure)
|
|---|
| 1422 | c jscl - how J scales
|
|---|
| 1423 | c 0: J = Jh
|
|---|
| 1424 | c 1: J = (L/2) Jh
|
|---|
| 1425 | c
|
|---|
| 1426 | c result is inexact because:
|
|---|
| 1427 | c neighbor's boundary condition at far end unknown
|
|---|
| 1428 | c length of neighbor's neighbor unknown
|
|---|
| 1429 | c (these contribs should be small for large N and
|
|---|
| 1430 | c elements of nearly equal size)
|
|---|
| 1431 | c
|
|---|
| 1432 | include 'SIZE'
|
|---|
| 1433 | real g(0:lx1-1,0:lx1-1)
|
|---|
| 1434 | real bh(0:lx1-1),jgl(1:lx2,0:lx1-1)
|
|---|
| 1435 | real ll,lm,lr
|
|---|
| 1436 | integer b0,b1
|
|---|
| 1437 | logical l,r
|
|---|
| 1438 | integer jscl
|
|---|
| 1439 | c
|
|---|
| 1440 | real bl(0:lx1-1),bm(0:lx1-1),br(0:lx1-1)
|
|---|
| 1441 | real gl,gm,gr,gll,glm,gmm,gmr,grr
|
|---|
| 1442 | real fac
|
|---|
| 1443 | integer n
|
|---|
| 1444 | n=lx1-1
|
|---|
| 1445 | c
|
|---|
| 1446 | c
|
|---|
| 1447 | c compute the scale factors for J
|
|---|
| 1448 | if (jscl.eq.0) then
|
|---|
| 1449 | gl=1.
|
|---|
| 1450 | gm=1.
|
|---|
| 1451 | gr=1.
|
|---|
| 1452 | elseif (jscl.eq.1) then
|
|---|
| 1453 | gl=0.5*ll
|
|---|
| 1454 | gm=0.5*lm
|
|---|
| 1455 | gr=0.5*lr
|
|---|
| 1456 | endif
|
|---|
| 1457 | gll = gl*gl
|
|---|
| 1458 | glm = gl*gm
|
|---|
| 1459 | gmm = gm*gm
|
|---|
| 1460 | gmr = gm*gr
|
|---|
| 1461 | grr = gr*gr
|
|---|
| 1462 | c
|
|---|
| 1463 | c compute the summed inverse mass matrices for
|
|---|
| 1464 | c the middle, left, and right elements
|
|---|
| 1465 | do i=1,n-1
|
|---|
| 1466 | bm(i)=2. /(lm*bh(i))
|
|---|
| 1467 | enddo
|
|---|
| 1468 | if (b0.eq.0) then
|
|---|
| 1469 | bm(0)=0.5*lm*bh(0)
|
|---|
| 1470 | if(l) bm(0)=bm(0)+0.5*ll*bh(n)
|
|---|
| 1471 | bm(0)=1. /bm(0)
|
|---|
| 1472 | endif
|
|---|
| 1473 | if (b1.eq.n) then
|
|---|
| 1474 | bm(n)=0.5*lm*bh(n)
|
|---|
| 1475 | if(r) bm(n)=bm(n)+0.5*lr*bh(0)
|
|---|
| 1476 | bm(n)=1. /bm(n)
|
|---|
| 1477 | endif
|
|---|
| 1478 | c note that in computing bl for the left element,
|
|---|
| 1479 | c bl(0) is missing the contribution from its left neighbor
|
|---|
| 1480 | if (l) then
|
|---|
| 1481 | do i=0,n-1
|
|---|
| 1482 | bl(i)=2. /(ll*bh(i))
|
|---|
| 1483 | enddo
|
|---|
| 1484 | bl(n)=bm(0)
|
|---|
| 1485 | endif
|
|---|
| 1486 | c note that in computing br for the right element,
|
|---|
| 1487 | c br(n) is missing the contribution from its right neighbor
|
|---|
| 1488 | if (r) then
|
|---|
| 1489 | do i=1,n
|
|---|
| 1490 | br(i)=2. /(lr*bh(i))
|
|---|
| 1491 | enddo
|
|---|
| 1492 | br(0)=bm(n)
|
|---|
| 1493 | endif
|
|---|
| 1494 | c
|
|---|
| 1495 | call rzero(g,(n+1)*(n+1))
|
|---|
| 1496 | do j=1,n-1
|
|---|
| 1497 | do i=1,n-1
|
|---|
| 1498 | do k=b0,b1
|
|---|
| 1499 | g(i,j) = g(i,j) + gmm*jgl(i,k)*bm(k)*jgl(j,k)
|
|---|
| 1500 | enddo
|
|---|
| 1501 | enddo
|
|---|
| 1502 | enddo
|
|---|
| 1503 | c
|
|---|
| 1504 | if (l) then
|
|---|
| 1505 | do i=1,n-1
|
|---|
| 1506 | g(i,0) = glm*jgl(i,0)*bm(0)*jgl(n-1,n)
|
|---|
| 1507 | g(0,i) = g(i,0)
|
|---|
| 1508 | enddo
|
|---|
| 1509 | c the following is inexact
|
|---|
| 1510 | c the neighbors bc's are ignored, and the contribution
|
|---|
| 1511 | c from the neighbor's neighbor is left out
|
|---|
| 1512 | c that is, bl(0) could be off as noted above
|
|---|
| 1513 | c or maybe i should go from 1 to n
|
|---|
| 1514 | do i=0,n
|
|---|
| 1515 | g(0,0) = g(0,0) + gll*jgl(n-1,i)*bl(i)*jgl(n-1,i)
|
|---|
| 1516 | enddo
|
|---|
| 1517 | else
|
|---|
| 1518 | g(0,0)=1.
|
|---|
| 1519 | endif
|
|---|
| 1520 | c
|
|---|
| 1521 | if (r) then
|
|---|
| 1522 | do i=1,n-1
|
|---|
| 1523 | g(i,n) = gmr*jgl(i,n)*bm(n)*jgl(1,0)
|
|---|
| 1524 | g(n,i) = g(i,n)
|
|---|
| 1525 | enddo
|
|---|
| 1526 | c the following is inexact
|
|---|
| 1527 | c the neighbors bc's are ignored, and the contribution
|
|---|
| 1528 | c from the neighbor's neighbor is left out
|
|---|
| 1529 | c that is, br(n) could be off as noted above
|
|---|
| 1530 | c or maybe i should go from 0 to n-1
|
|---|
| 1531 | do i=0,n
|
|---|
| 1532 | g(n,n) = g(n,n) + grr*jgl(1,i)*br(i)*jgl(1,i)
|
|---|
| 1533 | enddo
|
|---|
| 1534 | else
|
|---|
| 1535 | g(n,n)=1.
|
|---|
| 1536 | endif
|
|---|
| 1537 | return
|
|---|
| 1538 | end
|
|---|
| 1539 | c-----------------------------------------------------------------------
|
|---|
| 1540 | subroutine swap_lengths
|
|---|
| 1541 |
|
|---|
| 1542 | include 'SIZE'
|
|---|
| 1543 | include 'INPUT'
|
|---|
| 1544 | include 'GEOM'
|
|---|
| 1545 | include 'WZ'
|
|---|
| 1546 | common /swaplengths/ l(lx1,ly1,lz1,lelv)
|
|---|
| 1547 | common /ctmpf/ lr(2*lx1+4),ls(2*lx1+4),lt(2*lx1+4)
|
|---|
| 1548 | $ , llr(lelt),lls(lelt),llt(lelt)
|
|---|
| 1549 | $ , lmr(lelt),lms(lelt),lmt(lelt)
|
|---|
| 1550 | $ , lrr(lelt),lrs(lelt),lrt(lelt)
|
|---|
| 1551 | real lr ,ls ,lt
|
|---|
| 1552 | real llr,lls,llt
|
|---|
| 1553 | real lmr,lms,lmt
|
|---|
| 1554 | real lrr,lrs,lrt
|
|---|
| 1555 |
|
|---|
| 1556 | real l,l2d
|
|---|
| 1557 | integer e
|
|---|
| 1558 |
|
|---|
| 1559 | n2 = lx1-1
|
|---|
| 1560 | nz0 = 1
|
|---|
| 1561 | nzn = 1
|
|---|
| 1562 | nx = lx1-2
|
|---|
| 1563 | if (if3d) then
|
|---|
| 1564 | nz0 = 0
|
|---|
| 1565 | nzn = n2
|
|---|
| 1566 | endif
|
|---|
| 1567 | call plane_space(lmr,lms,lmt,0,n2,wxm1,xm1,ym1,zm1,nx,n2,nz0,nzn)
|
|---|
| 1568 |
|
|---|
| 1569 | n=n2+1
|
|---|
| 1570 | if (if3d) then
|
|---|
| 1571 | do e=1,nelv
|
|---|
| 1572 | do j=2,n2
|
|---|
| 1573 | do k=2,n2
|
|---|
| 1574 | l(1,k,j,e) = lmr(e)
|
|---|
| 1575 | l(n,k,j,e) = lmr(e)
|
|---|
| 1576 | l(k,1,j,e) = lms(e)
|
|---|
| 1577 | l(k,n,j,e) = lms(e)
|
|---|
| 1578 | l(k,j,1,e) = lmt(e)
|
|---|
| 1579 | l(k,j,n,e) = lmt(e)
|
|---|
| 1580 | enddo
|
|---|
| 1581 | enddo
|
|---|
| 1582 | enddo
|
|---|
| 1583 | call dssum(l,n,n,n)
|
|---|
| 1584 | do e=1,nelv
|
|---|
| 1585 | llr(e) = l(1,2,2,e)-lmr(e)
|
|---|
| 1586 | lrr(e) = l(n,2,2,e)-lmr(e)
|
|---|
| 1587 | lls(e) = l(2,1,2,e)-lms(e)
|
|---|
| 1588 | lrs(e) = l(2,n,2,e)-lms(e)
|
|---|
| 1589 | llt(e) = l(2,2,1,e)-lmt(e)
|
|---|
| 1590 | lrt(e) = l(2,2,n,e)-lmt(e)
|
|---|
| 1591 | enddo
|
|---|
| 1592 | else
|
|---|
| 1593 | do e=1,nelv
|
|---|
| 1594 | do j=2,n2
|
|---|
| 1595 | l(1,j,1,e) = lmr(e)
|
|---|
| 1596 | l(n,j,1,e) = lmr(e)
|
|---|
| 1597 | l(j,1,1,e) = lms(e)
|
|---|
| 1598 | l(j,n,1,e) = lms(e)
|
|---|
| 1599 | c call outmat(l(1,1,1,e),n,n,' L ',e)
|
|---|
| 1600 | enddo
|
|---|
| 1601 | enddo
|
|---|
| 1602 | c call outmat(l(1,1,1,25),n,n,' L ',25)
|
|---|
| 1603 | call dssum(l,n,n,1)
|
|---|
| 1604 | c call outmat(l(1,1,1,25),n,n,' L ',25)
|
|---|
| 1605 | do e=1,nelv
|
|---|
| 1606 | c call outmat(l(1,1,1,e),n,n,' L ',e)
|
|---|
| 1607 | llr(e) = l(1,2,1,e)-lmr(e)
|
|---|
| 1608 | lrr(e) = l(n,2,1,e)-lmr(e)
|
|---|
| 1609 | lls(e) = l(2,1,1,e)-lms(e)
|
|---|
| 1610 | lrs(e) = l(2,n,1,e)-lms(e)
|
|---|
| 1611 | enddo
|
|---|
| 1612 | endif
|
|---|
| 1613 | return
|
|---|
| 1614 | end
|
|---|
| 1615 | c----------------------------------------------------------------------
|
|---|
| 1616 | subroutine row_zero(a,m,n,e)
|
|---|
| 1617 | integer m,n,e
|
|---|
| 1618 | real a(m,n)
|
|---|
| 1619 | do j=1,n
|
|---|
| 1620 | a(e,j)=0.
|
|---|
| 1621 | enddo
|
|---|
| 1622 | return
|
|---|
| 1623 | end
|
|---|
| 1624 | c-----------------------------------------------------------------------
|
|---|
| 1625 | subroutine mhd_bc_dn(ibc,face,e)
|
|---|
| 1626 | integer face,e
|
|---|
| 1627 | c
|
|---|
| 1628 | c sets Neumann BC on pressure (ibc=2) for face and e(lement) except
|
|---|
| 1629 | c when ifield normal component has (homogeneous) Neumann
|
|---|
| 1630 | c boundary condition setting ibc=1 (i.e. Direchlet BC on pressure)
|
|---|
| 1631 | c
|
|---|
| 1632 | c Note: 'SYM' on a plane with r,s,t-normal is 'dnn','ndn','nnd'? bsym?
|
|---|
| 1633 | c
|
|---|
| 1634 | include 'SIZE'
|
|---|
| 1635 | include 'TOPOL'
|
|---|
| 1636 | include 'INPUT'
|
|---|
| 1637 | include 'TSTEP'
|
|---|
| 1638 |
|
|---|
| 1639 | ied = eface(face) ! symmetric -> preprocessor notation
|
|---|
| 1640 | nfc = face+1
|
|---|
| 1641 | nfc = nfc/2 ! = 1,2,3 for face 1 & 2,3 & 4,5 & 6
|
|---|
| 1642 |
|
|---|
| 1643 | if (indx1(cbc(ied,e,ifield),'d',1).gt.0) ibc=2
|
|---|
| 1644 |
|
|---|
| 1645 | if (indx1(cbc(ied,e,ifield),'n',1).gt.nfc) ibc=1 ! 'n' for V_n
|
|---|
| 1646 |
|
|---|
| 1647 | return
|
|---|
| 1648 | end
|
|---|
| 1649 | c-----------------------------------------------------------------------
|
|---|