Actual source code: ex54f.F90

  1: !
  2: !   Description: Solve Ax=b.  A comes from an anisotropic 2D thermal problem with Q1 FEM on domain (-1,1)^2.
  3: !       Material conductivity given by tensor:
  4: !
  5: !       D = | 1 0       |
  6: !           | 0 epsilon |
  7: !
  8: !    rotated by angle 'theta' (-theta <90> in degrees) with anisotropic parameter 'epsilon' (-epsilon <0.0>).
  9: !    Blob right hand side centered at C (-blob_center C(1),C(2) <0,0>)
 10: !    Dirichlet BCs on y=-1 face.
 11: !
 12: !    -out_matlab will generate binary files for A,x,b and a ex54f.m file that reads them and plots them in matlab.
 13: !
 14: !    User can change anisotropic shape with function ex54_psi().  Negative theta will switch to a circular anisotropy.
 15: !

 17: ! -----------------------------------------------------------------------
 18:       program main
 19: #include <petsc/finclude/petscksp.h>
 20:       use petscksp
 21:       implicit none

 23:       Vec              xvec,bvec,uvec
 24:       Mat              Amat
 25:       KSP              ksp
 26:       PetscErrorCode   ierr
 27:       PetscViewer viewer
 28:       PetscInt qj,qi,ne,M,Istart,Iend,geq
 29:       PetscInt ki,kj,nel,ll,j1,i1,ndf,f4
 30:       PetscInt f2,f9,f6,one
 31:       PetscInt :: idx(4)
 32:       PetscBool  flg,out_matlab
 33:       PetscMPIInt size,rank
 34:       PetscScalar::ss(4,4),val
 35:       PetscReal::shp(3,9),sg(3,9)
 36:       PetscReal::thk,a1,a2
 37:       PetscReal, external :: ex54_psi
 38:       PetscReal::theta,eps,h,x,y,xsj
 39:       PetscReal::coord(2,4),dd(2,2),ev(3),blb(2)

 41:       common /ex54_theta/ theta
 42: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 43: !                 Beginning of program
 44: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 45:       call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
 46:       if (ierr .ne. 0) then
 47:         print*,'Unable to initialize PETSc'
 48:         stop
 49:       endif
 50:       call MPI_Comm_size(PETSC_COMM_WORLD,size,ierr)
 51:       call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)
 52:       one = 1
 53: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 54: !                 set parameters
 55: !     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 56:       f4 = 4
 57:       f2 = 2
 58:       f9 = 9
 59:       f6 = 6
 60:       ne = 9
 61:       call PetscOptionsGetInt(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-ne',ne,flg,ierr)
 62:       h = 2.0/real(ne)
 63:       M = (ne+1)*(ne+1)
 64:       theta = 90.0
 65: !     theta is input in degrees
 66:       call PetscOptionsGetReal(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-theta',theta,flg,ierr)
 67:       theta = theta / 57.2957795
 68:       eps = 1.0
 69:       call PetscOptionsGetReal(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-epsilon',eps,flg,ierr)
 70:       ki = 2
 71:       call PetscOptionsGetRealArray(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-blob_center',blb,ki,flg,ierr)
 72:       if (.not. flg) then
 73:          blb(1) = 0.0
 74:          blb(2) = 0.0
 75:       else if (ki .ne. 2) then
 76:          print *, 'error: ', ki,' arguments read for -blob_center.  Needs to be two.'
 77:       endif
 78:       call PetscOptionsGetBool(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-out_matlab',out_matlab,flg,ierr)
 79:       if (.not.flg) out_matlab = PETSC_FALSE;

 81:       ev(1) = 1.0
 82:       ev(2) = eps*ev(1)
 83: !     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 84: !     Compute the matrix and right-hand-side vector that define
 85: !     the linear system, Ax = b.
 86: !     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 87: !  Create matrix.  When using MatCreate(), the matrix format can
 88: !  be specified at runtime.
 89:       call MatCreate(PETSC_COMM_WORLD,Amat,ierr)
 90:       call MatSetSizes( Amat,PETSC_DECIDE, PETSC_DECIDE, M, M, ierr)
 91:       call MatSetType( Amat, MATAIJ, ierr)
 92:       call MatSetOption(Amat,MAT_SPD,PETSC_TRUE,ierr)
 93:       if (size == 1) then
 94:          call MatSetType( Amat, MATAIJ, ierr)
 95:       else
 96:          call MatSetType( Amat, MATMPIAIJ, ierr)
 97:       endif
 98:       call MatMPIAIJSetPreallocation(Amat,f9,PETSC_NULL_INTEGER,f6,PETSC_NULL_INTEGER, ierr)
 99:       call MatSetFromOptions( Amat, ierr)
100:       call MatSetUp( Amat, ierr)
101:       call MatGetOwnershipRange( Amat, Istart, Iend, ierr)
102: !  Create vectors.  Note that we form 1 vector from scratch and
103: !  then duplicate as needed.
104:       call MatCreateVecs( Amat, PETSC_NULL_VEC, xvec, ierr)
105:       call VecSetFromOptions( xvec, ierr)
106:       call VecDuplicate( xvec, bvec, ierr)
107:       call VecDuplicate( xvec, uvec, ierr)
108: !  Assemble matrix.
109: !   - Note that MatSetValues() uses 0-based row and column numbers
110: !     in Fortran as well as in C (as set here in the array "col").
111:       thk = 1.0              ! thickness
112:       nel = 4                   ! nodes per element (quad)
113:       ndf = 1
114:       call int2d(f2,sg)
115:       do geq=Istart,Iend-1,1
116:          qj = geq/(ne+1); qi = mod(geq,(ne+1))
117:          x = h*qi - 1.0; y = h*qj - 1.0 ! lower left corner (-1,-1)
118:          if (qi < ne .and. qj < ne) then
119:             coord(1,1) = x;   coord(2,1) = y
120:             coord(1,2) = x+h; coord(2,2) = y
121:             coord(1,3) = x+h; coord(2,3) = y+h
122:             coord(1,4) = x;   coord(2,4) = y+h
123: ! form stiff
124:             ss = 0.0
125:             do ll = 1,4
126:                call shp2dquad(sg(1,ll),sg(2,ll),coord,shp,xsj,f2)
127:                xsj = xsj*sg(3,ll)*thk
128:                call thfx2d(ev,coord,shp,dd,f2,f2,f4,ex54_psi)
129:                j1 = 1
130:                do kj = 1,nel
131:                   a1 = (dd(1,1)*shp(1,kj) + dd(1,2)*shp(2,kj))*xsj
132:                   a2 = (dd(2,1)*shp(1,kj) + dd(2,2)*shp(2,kj))*xsj
133: !     Compute residual
134: !                  p(j1) = p(j1) - a1*gradt(1) - a2*gradt(2)
135: !     Compute tangent
136:                   i1 = 1
137:                   do ki = 1,nel
138:                      ss(i1,j1) = ss(i1,j1) + a1*shp(1,ki) + a2*shp(2,ki)
139:                      i1 = i1 + ndf
140:                   end do
141:                   j1 = j1 + ndf
142:                end do
143:             enddo

145:             idx(1) = geq; idx(2) = geq+1; idx(3) = geq+(ne+1)+1
146:             idx(4) = geq+(ne+1)
147:             if (qj > 0) then
148:                call MatSetValues(Amat,f4,idx,f4,idx,ss,ADD_VALUES,ierr)
149:             else                !     a BC
150:                do ki=1,4,1
151:                   do kj=1,4,1
152:                      if (ki<3 .or. kj<3) then
153:                         if (ki==kj) then
154:                            ss(ki,kj) = .1*ss(ki,kj)
155:                         else
156:                            ss(ki,kj) = 0.0
157:                         endif
158:                      endif
159:                   enddo
160:                enddo
161:                call MatSetValues(Amat,f4,idx,f4,idx,ss,ADD_VALUES,ierr)
162:             endif               ! BC
163:          endif                  ! add element
164:          if (qj > 0) then      ! set rhs
165:             val = h*h*exp(-100*((x+h/2)-blb(1))**2)*exp(-100*((y+h/2)-blb(2))**2)
166:             call VecSetValues(bvec,one,geq,val,INSERT_VALUES,ierr)
167:          endif
168:       enddo
169:       call MatAssemblyBegin(Amat,MAT_FINAL_ASSEMBLY,ierr)
170:       call MatAssemblyEnd(Amat,MAT_FINAL_ASSEMBLY,ierr)
171:       call VecAssemblyBegin(bvec,ierr)
172:       call VecAssemblyEnd(bvec,ierr)

174: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
175: !          Create the linear solver and set various options
176: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

178: !  Create linear solver context

180:       call KSPCreate(PETSC_COMM_WORLD,ksp,ierr)

182: !  Set operators. Here the matrix that defines the linear system
183: !  also serves as the preconditioning matrix.

185:       call KSPSetOperators(ksp,Amat,Amat,ierr)

187: !  Set runtime options, e.g.,
188: !      -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
189: !  These options will override those specified above as long as
190: !  KSPSetFromOptions() is called _after_ any other customization
191: !  routines.

193:       call KSPSetFromOptions(ksp,ierr)

195: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
196: !                      Solve the linear system
197: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

199:       call KSPSolve(ksp,bvec,xvec,ierr)
200:       CHKERRA(ierr)

202: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
203: !                      output
204: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
205:       if (out_matlab) then
206:          call PetscViewerBinaryOpen(PETSC_COMM_WORLD,'Amat',FILE_MODE_WRITE,viewer,ierr)
207:          call MatView(Amat,viewer,ierr)
208:          call PetscViewerDestroy(viewer,ierr)

210:          call PetscViewerBinaryOpen(PETSC_COMM_WORLD,'Bvec',FILE_MODE_WRITE,viewer,ierr)
211:          call VecView(bvec,viewer,ierr)
212:          call PetscViewerDestroy(viewer,ierr)

214:          call PetscViewerBinaryOpen(PETSC_COMM_WORLD,'Xvec',FILE_MODE_WRITE,viewer,ierr)
215:          call VecView(xvec,viewer,ierr)
216:          call PetscViewerDestroy(viewer,ierr)

218:          call MatMult(Amat,xvec,uvec,ierr)
219:          val = -1.0
220:          call VecAXPY(uvec,val,bvec,ierr)
221:          call PetscViewerBinaryOpen(PETSC_COMM_WORLD,'Rvec',FILE_MODE_WRITE,viewer,ierr)
222:          call VecView(uvec,viewer,ierr)
223:          call PetscViewerDestroy(viewer,ierr)

225:          if (rank == 0) then
226:             open(1,file='ex54f.m', FORM='formatted')
227:             write (1,*) 'A = PetscBinaryRead(''Amat'');'
228:             write (1,*) '[m n] = size(A);'
229:             write (1,*) 'mm = sqrt(m);'
230:             write (1,*) 'b = PetscBinaryRead(''Bvec'');'
231:             write (1,*) 'x = PetscBinaryRead(''Xvec'');'
232:             write (1,*) 'r = PetscBinaryRead(''Rvec'');'
233:             write (1,*) 'bb = reshape(b,mm,mm);'
234:             write (1,*) 'xx = reshape(x,mm,mm);'
235:             write (1,*) 'rr = reshape(r,mm,mm);'
236: !            write (1,*) 'imagesc(bb')'
237: !            write (1,*) 'title('RHS'),'
238:             write (1,*) 'figure,'
239:             write (1,*) 'imagesc(xx'')'
240:             write (1,2002) eps,theta*57.2957795
241:             write (1,*) 'figure,'
242:             write (1,*) 'imagesc(rr'')'
243:             write (1,*) 'title(''Residual''),'
244:             close(1)
245:          endif
246:       endif
247:  2002 format('title(''Solution: esp='',d9.3,'', theta='',g8.3,''),')
248: !  Free work space.  All PETSc objects should be destroyed when they
249: !  are no longer needed.

251:       call VecDestroy(xvec,ierr)
252:       call VecDestroy(bvec,ierr)
253:       call VecDestroy(uvec,ierr)
254:       call MatDestroy(Amat,ierr)
255:       call KSPDestroy(ksp,ierr)
256:       call PetscFinalize(ierr)

258:       end

260: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
261: !     thfx2d - compute material tensor
262: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
263: !     Compute thermal gradient and flux

265:       subroutine thfx2d(ev,xl,shp,dd,ndm,ndf,nel,dir)
266:       implicit  none

268:       PetscInt   ndm,ndf,nel,i
269:       PetscReal ev(2),xl(ndm,nel),shp(3,*),dir
270:       PetscReal xx,yy,psi,cs,sn,c2,s2,dd(2,2)

272:       xx       = 0.0
273:       yy       = 0.0
274:       do i = 1,nel
275:         xx       = xx       + shp(3,i)*xl(1,i)
276:         yy       = yy       + shp(3,i)*xl(2,i)
277:       end do
278:       psi = dir(xx,yy)
279: !     Compute thermal flux
280:       cs  = cos(psi)
281:       sn  = sin(psi)
282:       c2  = cs*cs
283:       s2  = sn*sn
284:       cs  = cs*sn

286:       dd(1,1) = c2*ev(1) + s2*ev(2)
287:       dd(2,2) = s2*ev(1) + c2*ev(2)
288:       dd(1,2) = cs*(ev(1) - ev(2))
289:       dd(2,1) = dd(1,2)

291: !      flux(1) = -dd(1,1)*gradt(1) - dd(1,2)*gradt(2)
292: !      flux(2) = -dd(2,1)*gradt(1) - dd(2,2)*gradt(2)

294:       end

296: !     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
297: !     shp2dquad - shape functions - compute derivatives w/r natural coords.
298: !     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
299:        subroutine shp2dquad(s,t,xl,shp,xsj,ndm)
300: !-----[--.----+----.----+----.-----------------------------------------]
301: !      Purpose: Shape function routine for 4-node isoparametric quads
302: !
303: !      Inputs:
304: !         s,t       - Natural coordinates of point
305: !         xl(ndm,*) - Nodal coordinates for element
306: !         ndm       - Spatial dimension of mesh

308: !      Outputs:
309: !         shp(3,*)  - Shape functions and derivatives at point
310: !                     shp(1,i) = dN_i/dx  or dN_i/dxi_1
311: !                     shp(2,i) = dN_i/dy  or dN_i/dxi_2
312: !                     shp(3,i) = N_i
313: !         xsj       - Jacobian determinant at point
314: !-----[--.----+----.----+----.-----------------------------------------]
315:       implicit  none
316:       PetscInt  ndm
317:       PetscReal xo,xs,xt, yo,ys,yt, xsm,xsp,xtm
318:       PetscReal xtp, ysm,ysp,ytm,ytp
319:       PetscReal s,t, xsj,xsj1, sh,th,sp,tp,sm
320:       PetscReal tm, xl(ndm,4),shp(3,4)

322: !     Set up interpolations

324:       sh = 0.5*s
325:       th = 0.5*t
326:       sp = 0.5 + sh
327:       tp = 0.5 + th
328:       sm = 0.5 - sh
329:       tm = 0.5 - th
330:       shp(3,1) =   sm*tm
331:       shp(3,2) =   sp*tm
332:       shp(3,3) =   sp*tp
333:       shp(3,4) =   sm*tp

335: !     Set up natural coordinate functions (times 4)

337:       xo =  xl(1,1)-xl(1,2)+xl(1,3)-xl(1,4)
338:       xs = -xl(1,1)+xl(1,2)+xl(1,3)-xl(1,4) + xo*t
339:       xt = -xl(1,1)-xl(1,2)+xl(1,3)+xl(1,4) + xo*s
340:       yo =  xl(2,1)-xl(2,2)+xl(2,3)-xl(2,4)
341:       ys = -xl(2,1)+xl(2,2)+xl(2,3)-xl(2,4) + yo*t
342:       yt = -xl(2,1)-xl(2,2)+xl(2,3)+xl(2,4) + yo*s

344: !     Compute jacobian (times 16)

346:       xsj1 = xs*yt - xt*ys

348: !     Divide jacobian by 16 (multiply by .0625)

350:       xsj = 0.0625*xsj1
351:       if (xsj1.eq.0.0) then
352:          xsj1 = 1.0
353:       else
354:          xsj1 = 1.0/xsj1
355:       endif

357: !     Divide functions by jacobian

359:       xs  = (xs+xs)*xsj1
360:       xt  = (xt+xt)*xsj1
361:       ys  = (ys+ys)*xsj1
362:       yt  = (yt+yt)*xsj1

364: !     Multiply by interpolations

366:       ytm =  yt*tm
367:       ysm =  ys*sm
368:       ytp =  yt*tp
369:       ysp =  ys*sp
370:       xtm =  xt*tm
371:       xsm =  xs*sm
372:       xtp =  xt*tp
373:       xsp =  xs*sp

375: !     Compute shape functions

377:       shp(1,1) = - ytm+ysm
378:       shp(1,2) =   ytm+ysp
379:       shp(1,3) =   ytp-ysp
380:       shp(1,4) = - ytp-ysm
381:       shp(2,1) =   xtm-xsm
382:       shp(2,2) = - xtm-xsp
383:       shp(2,3) = - xtp+xsp
384:       shp(2,4) =   xtp+xsm

386:       end

388: !     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
389: !     int2d
390: !     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
391:       subroutine int2d(l,sg)
392: !-----[--.----+----.----+----.-----------------------------------------]
393: !     Purpose: Form Gauss points and weights for two dimensions

395: !     Inputs:
396: !     l       - Number of points/direction

398: !     Outputs:
399: !     sg(3,*) - Array of points and weights
400: !-----[--.----+----.----+----.-----------------------------------------]
401:       implicit  none
402:       PetscInt   l,i,lr(9),lz(9)
403:       PetscReal    g,third,sg(3,*)
404:       data      lr/-1,1,1,-1,0,1,0,-1,0/,lz/-1,-1,1,1,-1,0,1,0,0/
405:       data      third / 0.3333333333333333 /

407: !     2x2 integration
408:       g = sqrt(third)
409:       do i = 1,4
410:          sg(1,i) = g*lr(i)
411:          sg(2,i) = g*lz(i)
412:          sg(3,i) = 1.0
413:       end do

415:       end

417: !     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
418: !     ex54_psi - anusotropic material direction
419: !     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
420:       PetscReal function ex54_psi(x,y)
421:       implicit  none
422:       PetscReal x,y,theta
423:       common /ex54_theta/ theta
424:       ex54_psi = theta
425:       if (theta < 0.) then     ! circular
426:          if (y==0) then
427:             ex54_psi = 2.0*atan(1.0)
428:          else
429:             ex54_psi = atan(-x/y)
430:          endif
431:       endif
432:       end

434: !
435: !/*TEST
436: !
437: !   build:
438: !
439: !   test:
440: !      nsize: 4
441: !      args: -ne 39 -theta 30.0 -epsilon 1.e-1 -blob_center 0.,0. -ksp_type cg -pc_type gamg -pc_gamg_type agg -pc_gamg_agg_nsmooths 1 -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 -mat_coarsen_type hem -pc_gamg_square_graph 0 -ksp_monitor_short -pc_gamg_esteig_ksp_max_it 5
442: !      requires: !single
443: !
444: !TEST*/