Actual source code: ex54f.F
petsc-3.3-p7 2013-05-11
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 <finclude/petscsys.h>
23: #include <finclude/petscvec.h>
24: #include <finclude/petscmat.h>
25: #include <finclude/petscksp.h>
26: #include <finclude/petscpc.h>
27: #include <finclude/petscviewer.h>
28: #include <finclude/petscviewer.h90>
29: Vec xvec,bvec,uvec
30: Mat Amat
31: KSP ksp
32: PC pc
33: PetscErrorCode ierr
34: PetscViewer viewer
35: PetscInt qj,qi,ne,M,Istart,Iend,geq,ix
36: PetscInt ki,kj,lint,nel,ll,j1,i1,ndf
37: PetscInt :: idx(4)
38: PetscBool flg,out_matlab
39: PetscMPIInt npe,mype
40: PetscScalar::ss(4,4),res(4),val
41: PetscReal::shp(3,9),sg(3,9),flux(2)
42: PetscReal::thk,a1,a2
43: PetscReal, external :: ex54_psi
44: PetscReal::norm,tol,theta,eps,h,x,y,xsj
45: PetscReal::coord(2,4),dd(2,2),ev(3),blb(2)
46:
47: common /ex54_theta/ theta
48: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
49: ! Beginning of program
50: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
51: call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
52: call MPI_Comm_size(PETSC_COMM_WORLD,npe,ierr)
53: call MPI_Comm_rank(PETSC_COMM_WORLD,mype,ierr)
54: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
55: ! set parameters
56: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
57: ne = 9
58: call PetscOptionsGetInt(PETSC_NULL_CHARACTER,'-ne',ne,flg,ierr)
59: h = 2.d0/ne
60: M = (ne+1)*(ne+1)
61: theta = 90.d0
62: ! theta is input in degrees
63: call PetscOptionsGetReal(PETSC_NULL_CHARACTER,'-theta',theta,flg,
64: $ ierr)
65: theta = theta / 57.2957795
66: eps = 1.d0
67: call PetscOptionsGetReal(PETSC_NULL_CHARACTER,'-epsilon',eps,flg,
68: $ ierr)
69: ki = 2
70: call PetscOptionsGetRealArray(PETSC_NULL_CHARACTER,'-blob_center',
71: $ blb,ki,flg,ierr)
72: if( .not. flg ) then
73: blb(1) = 0.d0
74: blb(2) = 0.d0
75: else if( ki .ne. 2 ) then
76: print *, 'error: ', ki,
77: $ ' arguments read for -blob_center. Needs to be two.'
78: endif
79: call PetscOptionsGetBool(PETSC_NULL_CHARACTER,'-out_matlab',
80: $ out_matlab,flg,ierr)
82: ev(1) = 1.d0
83: ev(2) = eps*ev(1)
84: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
85: ! Compute the matrix and right-hand-side vector that define
86: ! the linear system, Ax = b.
87: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
88: ! Create matrix. When using MatCreate(), the matrix format can
89: ! be specified at runtime.
90: call MatCreate(PETSC_COMM_WORLD,Amat,ierr)
91: call MatSetSizes( Amat,PETSC_DECIDE, PETSC_DECIDE, M, M, ierr )
92: if( npe == 1 ) then
93: call MatSetType( Amat, MATAIJ, ierr )
94: else
95: call MatSetType( Amat, MATMPIAIJ, ierr )
96: endif
97: call MatMPIAIJSetPreallocation(Amat,9,PETSC_NULL_INTEGER,6,
98: $ 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 MatGetVecs( Amat, PETSC_NULL_OBJECT, 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.d0 ! thickness
112: nel = 4 ! nodes per element (quad)
113: ndf = 1
114: call int2d(2,lint,sg)
115: ix = 0
116: do geq=Istart,Iend-1,1
117: qj = geq/(ne+1); qi = mod(geq,(ne+1))
118: x = h*qi - 1.d0; y = h*qj - 1.d0 ! lower left corner (-1,-1)
119: if( qi < ne .and. qj < ne ) then
120: coord(1,1) = x; coord(2,1) = y
121: coord(1,2) = x+h; coord(2,2) = y
122: coord(1,3) = x+h; coord(2,3) = y+h
123: coord(1,4) = x; coord(2,4) = y+h
124: ! form stiff
125: ss = 0.d0
126: do ll = 1,lint
127: call shp2dquad(sg(1,ll),sg(2,ll),coord,shp,xsj,2)
128: xsj = xsj*sg(3,ll)*thk
129: call thfx2d(ev,coord,shp,dd,2,2,4,ex54_psi)
130: j1 = 1
131: do kj = 1,nel
132: a1 = (dd(1,1)*shp(1,kj) + dd(1,2)*shp(2,kj))*xsj
133: a2 = (dd(2,1)*shp(1,kj) + dd(2,2)*shp(2,kj))*xsj
134: c Compute residual
135: c p(j1) = p(j1) - a1*gradt(1) - a2*gradt(2)
136: c Compute tangent
137: i1 = 1
138: do ki = 1,nel
139: ss(i1,j1) = ss(i1,j1) + a1*shp(1,ki) + a2*shp(2,ki)
140: i1 = i1 + ndf
141: end do
142: j1 = j1 + ndf
143: end do
144: enddo
146: idx(1) = geq; idx(2) = geq+1; idx(3) = geq+(ne+1)+1
147: idx(4) = geq+(ne+1)
148: if( qj > 0 ) then
149: call MatSetValues(Amat,4,idx,4,idx,ss,ADD_VALUES,ierr)
150: else ! a BC
151: do ki=1,4,1
152: do kj=1,4,1
153: if(ki<3 .or. kj<3 ) then
154: if( ki==kj ) then
155: ss(ki,kj) = .1d0*ss(ki,kj)
156: else
157: ss(ki,kj) = 0.d0
158: endif
159: endif
160: enddo
161: enddo
162: call MatSetValues(Amat,4,idx,4,idx,ss,ADD_VALUES,ierr)
163: endif ! BC
164: endif ! add element
165: if( qj > 0 ) then ! set rhs
167: val = h*h*exp(-1.d2*((x+h/2)-blb(1))**2)*
168: $ exp(-1.d2*((y+h/2)-blb(2))**2)
169: call VecSetValues(bvec,1,geq,val,INSERT_VALUES,ierr)
170: endif
171: enddo
172: call MatAssemblyBegin(Amat,MAT_FINAL_ASSEMBLY,ierr)
173: call MatAssemblyEnd(Amat,MAT_FINAL_ASSEMBLY,ierr)
174: call VecAssemblyBegin(bvec,ierr)
175: call VecAssemblyEnd(bvec,ierr)
177: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
178: ! Create the linear solver and set various options
179: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
181: ! Create linear solver context
183: call KSPCreate(PETSC_COMM_WORLD,ksp,ierr)
185: ! Set operators. Here the matrix that defines the linear system
186: ! also serves as the preconditioning matrix.
188: call KSPSetOperators(ksp,Amat,Amat,DIFFERENT_NONZERO_PATTERN,ierr)
190: ! Set runtime options, e.g.,
191: ! -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
192: ! These options will override those specified above as long as
193: ! KSPSetFromOptions() is called _after_ any other customization
194: ! routines.
196: call KSPSetFromOptions(ksp,ierr)
198: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
199: ! Solve the linear system
200: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
202: call KSPSolve(ksp,bvec,xvec,ierr)
205: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
206: ! output
207: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
208: if( out_matlab ) then
209: call PetscViewerBinaryOpen(PETSC_COMM_WORLD,'Amat',
210: $ FILE_MODE_WRITE,viewer,ierr)
211: call MatView(Amat,viewer,ierr)
212: call PetscViewerDestroy(viewer,ierr)
214: call PetscViewerBinaryOpen(PETSC_COMM_WORLD,'Bvec',
215: $ FILE_MODE_WRITE,viewer,ierr)
216: call VecView(bvec,viewer,ierr)
217: call PetscViewerDestroy(viewer,ierr)
219: call PetscViewerBinaryOpen(PETSC_COMM_WORLD,'Xvec',
220: $ FILE_MODE_WRITE,viewer,ierr)
221: call VecView(xvec,viewer,ierr)
222: call PetscViewerDestroy(viewer,ierr)
224: call MatMult(Amat,xvec,uvec,ierr)
225: val = -1.d0
226: call VecAXPY(uvec,val,bvec,ierr)
227: call PetscViewerBinaryOpen(PETSC_COMM_WORLD,'Rvec',
228: $ FILE_MODE_WRITE,viewer,ierr)
229: call VecView(uvec,viewer,ierr)
230: call PetscViewerDestroy(viewer,ierr)
232: if( mype == 0 ) then
233: open(1,file="ex54f.m", FORM="formatted")
234: write (1,*) "A = PetscBinaryRead('Amat');"
235: write (1,*) "[m n] = size(A);"
236: write (1,*) "mm = sqrt(m);"
237: write (1,*) "b = PetscBinaryRead('Bvec');"
238: write (1,*) "x = PetscBinaryRead('Xvec');"
239: write (1,*) "r = PetscBinaryRead('Rvec');"
240: write (1,*) "bb = reshape(b,mm,mm);"
241: write (1,*) "xx = reshape(x,mm,mm);"
242: write (1,*) "rr = reshape(r,mm,mm);"
243: c write (1,*) "imagesc(bb')"
244: c write (1,*) "title('RHS'),"
245: write (1,*) "figure,"
246: write (1,*) "imagesc(xx')"
247: write (1,2002) eps,theta*57.2957795
248: write (1,*) "figure,"
249: write (1,*) "imagesc(rr')"
250: write (1,*) "title('Residual'),"
251: close(1)
252: endif
253: endif
254: 2002 format("title('Solution: esp=",d9.3,", theta=",g8.3,"'),")
255: ! Free work space. All PETSc objects should be destroyed when they
256: ! are no longer needed.
258: call VecDestroy(xvec,ierr)
259: call VecDestroy(bvec,ierr)
260: call VecDestroy(uvec,ierr)
261: call MatDestroy(Amat,ierr)
262: call KSPDestroy(ksp,ierr)
263: call PetscFinalize(ierr)
265: end
267: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
268: ! thfx2d - compute material tensor
269: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
271: subroutine thfx2d(ev,xl,shp,dd,ndm,ndf,nel,dir)
272:
273: c Compute thermal gradient and flux
274:
275: implicit none
277: integer ndm,ndf,nel,i
278: PetscReal ev(2),xl(ndm,nel),shp(3,*),dir
279: PetscReal xx,yy,psi,cs,sn,c2,s2,dd(2,2)
281: c temp = 0.0d0
282: c gradt(1) = 0.0d0
283: c gradt(2) = 0.0d0
284: xx = 0.0d0
285: yy = 0.0d0
286: do i = 1,nel
287: c gradt(1) = gradt(1) + shp(1,i)*ul(1,i)
288: c gradt(2) = gradt(2) + shp(2,i)*ul(1,i)
289: c temp = temp + shp(3,i)*ul(1,i)
290: xx = xx + shp(3,i)*xl(1,i)
291: yy = yy + shp(3,i)*xl(2,i)
292: end do
293: psi = dir(xx,yy)
294: c Compute thermal flux
295: cs = cos(psi)
296: sn = sin(psi)
297: c2 = cs*cs
298: s2 = sn*sn
299: cs = cs*sn
301: dd(1,1) = c2*ev(1) + s2*ev(2)
302: dd(2,2) = s2*ev(1) + c2*ev(2)
303: dd(1,2) = cs*(ev(1) - ev(2))
304: dd(2,1) = dd(1,2)
306: c flux(1) = -dd(1,1)*gradt(1) - dd(1,2)*gradt(2)
307: c flux(2) = -dd(2,1)*gradt(1) - dd(2,2)*gradt(2)
309: end
311: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
312: ! shp2dquad - shape functions - compute derivatives w/r natural coords.
313: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
314: subroutine shp2dquad(s,t,xl,shp,xsj,ndm)
315: c-----[--.----+----.----+----.-----------------------------------------]
316: c Purpose: Shape function routine for 4-node isoparametric quads
317: c
318: c Inputs:
319: c s,t - Natural coordinates of point
320: c xl(ndm,*) - Nodal coordinates for element
321: c ndm - Spatial dimension of mesh
323: c Outputs:
324: c shp(3,*) - Shape functions and derivatives at point
325: c shp(1,i) = dN_i/dx or dN_i/dxi_1
326: c shp(2,i) = dN_i/dy or dN_i/dxi_2
327: c shp(3,i) = N_i
328: c xsj - Jacobian determinant at point
329: c-----[--.----+----.----+----.-----------------------------------------]
330: implicit none
331: integer ndm
332: real*8 xo,xs,xt, yo,ys,yt, xsm,xsp,xtm,xtp, ysm,ysp,ytm,ytp
333: real*8 s,t, xsj,xsj1, sh,th,sp,tp,sm,tm, xl(ndm,4),shp(3,4)
335: c Set up interpolations
337: sh = 0.5d0*s
338: th = 0.5d0*t
339: sp = 0.5d0 + sh
340: tp = 0.5d0 + th
341: sm = 0.5d0 - sh
342: tm = 0.5d0 - th
343: shp(3,1) = sm*tm
344: shp(3,2) = sp*tm
345: shp(3,3) = sp*tp
346: shp(3,4) = sm*tp
348: c Set up natural coordinate functions (times 4)
350: xo = xl(1,1)-xl(1,2)+xl(1,3)-xl(1,4)
351: xs = -xl(1,1)+xl(1,2)+xl(1,3)-xl(1,4) + xo*t
352: xt = -xl(1,1)-xl(1,2)+xl(1,3)+xl(1,4) + xo*s
353: yo = xl(2,1)-xl(2,2)+xl(2,3)-xl(2,4)
354: ys = -xl(2,1)+xl(2,2)+xl(2,3)-xl(2,4) + yo*t
355: yt = -xl(2,1)-xl(2,2)+xl(2,3)+xl(2,4) + yo*s
357: c Compute jacobian (times 16)
359: xsj1 = xs*yt - xt*ys
361: c Divide jacobian by 16 (multiply by .0625)
363: xsj = 0.0625d0*xsj1
364: if(xsj1.eq.0.0d0) then
365: xsj1 = 1.0d0
366: else
367: xsj1 = 1.0d0/xsj1
368: endif
369:
370: c Divide functions by jacobian
371:
372: xs = (xs+xs)*xsj1
373: xt = (xt+xt)*xsj1
374: ys = (ys+ys)*xsj1
375: yt = (yt+yt)*xsj1
376:
377: c Multiply by interpolations
378:
379: ytm = yt*tm
380: ysm = ys*sm
381: ytp = yt*tp
382: ysp = ys*sp
383: xtm = xt*tm
384: xsm = xs*sm
385: xtp = xt*tp
386: xsp = xs*sp
387:
388: c Compute shape functions
389:
390: shp(1,1) = - ytm+ysm
391: shp(1,2) = ytm+ysp
392: shp(1,3) = ytp-ysp
393: shp(1,4) = - ytp-ysm
394: shp(2,1) = xtm-xsm
395: shp(2,2) = - xtm-xsp
396: shp(2,3) = - xtp+xsp
397: shp(2,4) = xtp+xsm
398:
399: end
401: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
402: ! int2d
403: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
404: subroutine int2d(l,lint,sg)
405: c-----[--.----+----.----+----.-----------------------------------------]
406: c Purpose: Form Gauss points and weights for two dimensions
407:
408: c Inputs:
409: c l - Number of points/direction
410:
411: c Outputs:
412: c lint - Total number of points
413: c sg(3,*) - Array of points and weights
414: c-----[--.----+----.----+----.-----------------------------------------]
415: implicit none
416: integer l,i,lint,lr(9),lz(9)
417: real*8 g,third,sg(3,*)
418: data lr/-1,1,1,-1,0,1,0,-1,0/,lz/-1,-1,1,1,-1,0,1,0,0/
419: data third / 0.3333333333333333d0 /
420:
421: c Set number of total points
422:
423: lint = l*l
424:
425: c 2x2 integration
426: g = sqrt(third)
427: do i = 1,4
428: sg(1,i) = g*lr(i)
429: sg(2,i) = g*lz(i)
430: sg(3,i) = 1.d0
431: end do
433: end
435: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
436: ! ex54_psi - anusotropic material direction
437: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
438: PetscReal function ex54_psi(x,y)
439: implicit none
440: PetscReal x,y,theta
441: common /ex54_theta/ theta
442: ex54_psi = theta
443: if( theta < 0. ) then ! circular
444: if(y==0) then
445: ex54_psi = 2.d0*atan(1.d0)
446: else
447: ex54_psi = atan(-x/y)
448: endif
449: endif
450: end