Actual source code: ex54f.F

petsc-master 2016-12-04
Report Typos and Errors
  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: !
 16: !/*T
 17: !   Concepts: KSP^solving a system of linear equations
 18: !T*/
 19: ! -----------------------------------------------------------------------
 20:       program main
 21:       implicit none
 22:  #include <petsc/finclude/petscsys.h>
 23:  #include <petsc/finclude/petscvec.h>
 24:  #include <petsc/finclude/petscmat.h>
 25:  #include <petsc/finclude/petscksp.h>
 26:  #include <petsc/finclude/petscpc.h>
 27:  #include <petsc/finclude/petscviewer.h>
 28: #include <petsc/finclude/petscviewer.h90>
 29:       Vec              xvec,bvec,uvec
 30:       Mat              Amat
 31:       KSP              ksp
 32:       PetscErrorCode   ierr
 33:       PetscViewer viewer
 34:       PetscInt qj,qi,ne,M,Istart,Iend,geq,ix
 35:       PetscInt ki,kj,lint,nel,ll,j1,i1,ndf,f4
 36:       PetscInt f2,f9,f6
 37:       PetscInt :: idx(4)
 38:       PetscBool  flg,out_matlab
 39:       PetscMPIInt size,rank
 40:       PetscScalar::ss(4,4),val
 41:       PetscReal::shp(3,9),sg(3,9)
 42:       PetscReal::thk,a1,a2
 43:       PetscReal, external :: ex54_psi
 44:       PetscReal::theta,eps,h,x,y,xsj
 45:       PetscReal::coord(2,4),dd(2,2),ev(3),blb(2)

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

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

157:             idx(1) = geq; idx(2) = geq+1; idx(3) = geq+(ne+1)+1
158:             idx(4) = geq+(ne+1)
159:             if ( qj > 0 ) then
160:                call MatSetValues(Amat,f4,idx,f4,idx,ss,ADD_VALUES,ierr)
161:             else                !     a BC
162:                do ki=1,4,1
163:                   do kj=1,4,1
164:                      if (ki<3 .or. kj<3 ) then
165:                         if ( ki==kj ) then
166:                            ss(ki,kj) = .1*ss(ki,kj)
167:                         else
168:                            ss(ki,kj) = 0.0
169:                         endif
170:                      endif
171:                   enddo
172:                enddo
173:                call MatSetValues(Amat,f4,idx,f4,idx,ss,ADD_VALUES,ierr)
174:             endif               ! BC
175:          endif                  ! add element
176:          if ( qj > 0 ) then      ! set rhs

178:             val = h*h*exp(-100.*((x+h/2)-blb(1))**2)*                            &
179:      &           exp(-100*((y+h/2)-blb(2))**2)
180:             call VecSetValues(bvec,1,geq,val,INSERT_VALUES,ierr)
181:          endif
182:       enddo
183:       call MatAssemblyBegin(Amat,MAT_FINAL_ASSEMBLY,ierr)
184:       call MatAssemblyEnd(Amat,MAT_FINAL_ASSEMBLY,ierr)
185:       call VecAssemblyBegin(bvec,ierr)
186:       call VecAssemblyEnd(bvec,ierr)

188: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
189: !          Create the linear solver and set various options
190: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

192: !  Create linear solver context

194:       call KSPCreate(PETSC_COMM_WORLD,ksp,ierr)

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

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

201: !  Set runtime options, e.g.,
202: !      -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
203: !  These options will override those specified above as long as
204: !  KSPSetFromOptions() is called _after_ any other customization
205: !  routines.

207:       call KSPSetFromOptions(ksp,ierr)

209: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
210: !                      Solve the linear system
211: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

213:       call KSPSolve(ksp,bvec,xvec,ierr)
214:       CHKERRQ(ierr)


217: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
218: !                      output
219: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
220:       if ( out_matlab ) then
221:          call PetscViewerBinaryOpen(PETSC_COMM_WORLD,'Amat',                &
222:      &        FILE_MODE_WRITE,viewer,ierr)
223:          call MatView(Amat,viewer,ierr)
224:          call PetscViewerDestroy(viewer,ierr)

226:          call PetscViewerBinaryOpen(PETSC_COMM_WORLD,'Bvec',                &
227:      &        FILE_MODE_WRITE,viewer,ierr)
228:          call VecView(bvec,viewer,ierr)
229:          call PetscViewerDestroy(viewer,ierr)

231:          call PetscViewerBinaryOpen(PETSC_COMM_WORLD,'Xvec',                &
232:      &        FILE_MODE_WRITE,viewer,ierr)
233:          call VecView(xvec,viewer,ierr)
234:          call PetscViewerDestroy(viewer,ierr)

236:          call MatMult(Amat,xvec,uvec,ierr)
237:          val = -1.0
238:          call VecAXPY(uvec,val,bvec,ierr)
239:          call PetscViewerBinaryOpen(PETSC_COMM_WORLD,'Rvec',                &
240:      &        FILE_MODE_WRITE,viewer,ierr)
241:          call VecView(uvec,viewer,ierr)
242:          call PetscViewerDestroy(viewer,ierr)

244:          if ( rank == 0 ) then
245:             open(1,file='ex54f.m', FORM='formatted')
246:             write (1,*) 'A = PetscBinaryRead(''Amat'');'
247:             write (1,*) '[m n] = size(A);'
248:             write (1,*) 'mm = sqrt(m);'
249:             write (1,*) 'b = PetscBinaryRead(''Bvec'');'
250:             write (1,*) 'x = PetscBinaryRead(''Xvec'');'
251:             write (1,*) 'r = PetscBinaryRead(''Rvec'');'
252:             write (1,*) 'bb = reshape(b,mm,mm);'
253:             write (1,*) 'xx = reshape(x,mm,mm);'
254:             write (1,*) 'rr = reshape(r,mm,mm);'
255: !            write (1,*) 'imagesc(bb')'
256: !            write (1,*) 'title('RHS'),'
257:             write (1,*) 'figure,'
258:             write (1,*) 'imagesc(xx'')'
259:             write (1,2002) eps,theta*57.2957795
260:             write (1,*) 'figure,'
261:             write (1,*) 'imagesc(rr'')'
262:             write (1,*) 'title(''Residual''),'
263:             close(1)
264:          endif
265:       endif
266:  2002 format('title(''Solution: esp='',d9.3,'', theta='',g8.3,''),')
267: !  Free work space.  All PETSc objects should be destroyed when they
268: !  are no longer needed.

270:       call VecDestroy(xvec,ierr)
271:       call VecDestroy(bvec,ierr)
272:       call VecDestroy(uvec,ierr)
273:       call MatDestroy(Amat,ierr)
274:       call KSPDestroy(ksp,ierr)
275:       call PetscFinalize(ierr)

277:       end

279: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
280: !     thfx2d - compute material tensor
281: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

283:       subroutine thfx2d(ev,xl,shp,dd,ndm,ndf,nel,dir)

285: !     Compute thermal gradient and flux

287:       implicit  none

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

293:       xx       = 0.0
294:       yy       = 0.0
295:       do i = 1,nel
296:         xx       = xx       + shp(3,i)*xl(1,i)
297:         yy       = yy       + shp(3,i)*xl(2,i)
298:       end do
299:       psi = dir(xx,yy)
300: !     Compute thermal flux
301:       cs  = cos(psi)
302:       sn  = sin(psi)
303:       c2  = cs*cs
304:       s2  = sn*sn
305:       cs  = cs*sn

307:       dd(1,1) = c2*ev(1) + s2*ev(2)
308:       dd(2,2) = s2*ev(1) + c2*ev(2)
309:       dd(1,2) = cs*(ev(1) - ev(2))
310:       dd(2,1) = dd(1,2)

312: !      flux(1) = -dd(1,1)*gradt(1) - dd(1,2)*gradt(2)
313: !      flux(2) = -dd(2,1)*gradt(1) - dd(2,2)*gradt(2)

315:       end

317: !     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
318: !     shp2dquad - shape functions - compute derivatives w/r natural coords.
319: !     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
320:        subroutine shp2dquad(s,t,xl,shp,xsj,ndm)
321: !-----[--.----+----.----+----.-----------------------------------------]
322: !      Purpose: Shape function routine for 4-node isoparametric quads
323: !
324: !      Inputs:
325: !         s,t       - Natural coordinates of point
326: !         xl(ndm,*) - Nodal coordinates for element
327: !         ndm       - Spatial dimension of mesh

329: !      Outputs:
330: !         shp(3,*)  - Shape functions and derivatives at point
331: !                     shp(1,i) = dN_i/dx  or dN_i/dxi_1
332: !                     shp(2,i) = dN_i/dy  or dN_i/dxi_2
333: !                     shp(3,i) = N_i
334: !         xsj       - Jacobian determinant at point
335: !-----[--.----+----.----+----.-----------------------------------------]
336:       implicit  none
337:       PetscInt  ndm
338:       PetscReal xo,xs,xt, yo,ys,yt, xsm,xsp,xtm
339:       PetscReal xtp, ysm,ysp,ytm,ytp
340:       PetscReal s,t, xsj,xsj1, sh,th,sp,tp,sm
341:       PetscReal tm, xl(ndm,4),shp(3,4)

343: !     Set up interpolations

345:       sh = 0.5*s
346:       th = 0.5*t
347:       sp = 0.5 + sh
348:       tp = 0.5 + th
349:       sm = 0.5 - sh
350:       tm = 0.5 - th
351:       shp(3,1) =   sm*tm
352:       shp(3,2) =   sp*tm
353:       shp(3,3) =   sp*tp
354:       shp(3,4) =   sm*tp

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

358:       xo =  xl(1,1)-xl(1,2)+xl(1,3)-xl(1,4)
359:       xs = -xl(1,1)+xl(1,2)+xl(1,3)-xl(1,4) + xo*t
360:       xt = -xl(1,1)-xl(1,2)+xl(1,3)+xl(1,4) + xo*s
361:       yo =  xl(2,1)-xl(2,2)+xl(2,3)-xl(2,4)
362:       ys = -xl(2,1)+xl(2,2)+xl(2,3)-xl(2,4) + yo*t
363:       yt = -xl(2,1)-xl(2,2)+xl(2,3)+xl(2,4) + yo*s

365: !     Compute jacobian (times 16)

367:       xsj1 = xs*yt - xt*ys

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

371:       xsj = 0.0625*xsj1
372:       if (xsj1.eq.0.0) then
373:          xsj1 = 1.0
374:       else
375:          xsj1 = 1.0/xsj1
376:       endif

378: !     Divide functions by jacobian

380:       xs  = (xs+xs)*xsj1
381:       xt  = (xt+xt)*xsj1
382:       ys  = (ys+ys)*xsj1
383:       yt  = (yt+yt)*xsj1

385: !     Multiply by interpolations

387:       ytm =  yt*tm
388:       ysm =  ys*sm
389:       ytp =  yt*tp
390:       ysp =  ys*sp
391:       xtm =  xt*tm
392:       xsm =  xs*sm
393:       xtp =  xt*tp
394:       xsp =  xs*sp

396: !     Compute shape functions

398:       shp(1,1) = - ytm+ysm
399:       shp(1,2) =   ytm+ysp
400:       shp(1,3) =   ytp-ysp
401:       shp(1,4) = - ytp-ysm
402:       shp(2,1) =   xtm-xsm
403:       shp(2,2) = - xtm-xsp
404:       shp(2,3) = - xtp+xsp
405:       shp(2,4) =   xtp+xsm

407:       end

409: !     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
410: !     int2d
411: !     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
412:       subroutine int2d(l,sg)
413: !-----[--.----+----.----+----.-----------------------------------------]
414: !     Purpose: Form Gauss points and weights for two dimensions

416: !     Inputs:
417: !     l       - Number of points/direction

419: !     Outputs:
420: !     lint    - Total number of points
421: !     sg(3,*) - Array of points and weights
422: !-----[--.----+----.----+----.-----------------------------------------]
423:       implicit  none
424:       PetscInt   l,i,lint,lr(9),lz(9)
425:       PetscReal    g,third,sg(3,*)
426:       data      lr/-1,1,1,-1,0,1,0,-1,0/,lz/-1,-1,1,1,-1,0,1,0,0/
427:       data      third / 0.3333333333333333 /

429: !     Set number of total points

431:       lint = l*l

433: !     2x2 integration
434:       g = sqrt(third)
435:       do i = 1,4
436:          sg(1,i) = g*lr(i)
437:          sg(2,i) = g*lz(i)
438:          sg(3,i) = 1.0
439:       end do

441:       end

443: !     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
444: !     ex54_psi - anusotropic material direction
445: !     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
446:       PetscReal function ex54_psi(x,y)
447:       implicit  none
448:       PetscReal x,y,theta
449:       common /ex54_theta/ theta
450:       ex54_psi = theta
451:       if ( theta < 0. ) then     ! circular
452:          if (y==0) then
453:             ex54_psi = 2.0*atan(1.0)
454:          else
455:             ex54_psi = atan(-x/y)
456:          endif
457:       endif
458:       end