PageRenderTime 27ms CodeModel.GetById 20ms RepoModel.GetById 0ms app.codeStats 0ms

/src/ksp/ksp/examples/tutorials/ex54f.F.html

https://bitbucket.org/alexei-matveev/petsc-debian-pkg
HTML | 465 lines | 422 code | 43 blank | 0 comment | 0 complexity | e351e478a7c1f3d51c67f352dc9c2daa MD5 | raw file
  1. <center><a href="ex54f.F">Actual source code: ex54f.F</a></center><br>
  2. <html>
  3. <head> <link rel="canonical" href="http://www.mcs.anl.gov/petsc/petsc-current/src/ksp/ksp/examples/tutorials/ex54f.F.html" />
  4. <title></title>
  5. <meta name="generator" content="c2html 0.9.5">
  6. <meta name="date" content="2012-08-29T17:07:05+00:00">
  7. </head>
  8. <body bgcolor="#FFFFFF">
  9. <div id="version" align=right><b>petsc-3.3-p3 2012-08-29</b></div>
  10. <pre width="80"><a name="line1"> 1: </a>!
  11. <a name="line2"> 2: </a>! Description: Solve Ax=b. A comes from an anisotropic 2D thermal problem with Q1 FEM on domain (-1,1)^2.
  12. <a name="line3"> 3: </a>! Material conductivity given by tensor:
  13. <a name="line4"> 4: </a>!
  14. <a name="line5"> 5: </a>! D = | 1 0 |
  15. <a name="line6"> 6: </a>! | 0 epsilon |
  16. <a name="line7"> 7: </a>!
  17. <a name="line8"> 8: </a>! rotated by angle 'theta' (-theta &lt;90&gt; in degrees) with anisotropic parameter 'epsilon' (-epsilon &lt;0.0&gt;).
  18. <a name="line9"> 9: </a>! Blob right hand side centered at C (-blob_center C(1),C(2) &lt;0,0&gt;)
  19. <a name="line10"> 10: </a>! Dirichlet BCs on y=-1 face.
  20. <a name="line11"> 11: </a>!
  21. <a name="line12"> 12: </a>! -out_matlab will generate binary files <font color="#4169E1">for</font> A,x,b and a ex54f.m file that reads them and plots them in matlab.
  22. <a name="line13"> 13: </a>!
  23. <a name="line14"> 14: </a>! User can change anisotropic shape with function ex54_psi(). Negative theta will <font color="#4169E1">switch</font> to a circular anisotropy.
  24. <a name="line15"> 15: </a>!
  25. <a name="line16"> 16: </a>!<font color="#B22222">/*T</font>
  26. <a name="line17"> 17: </a><font color="#B22222">! Concepts: <A href="../../../../../docs/manualpages/KSP/KSP.html#KSP">KSP</A>^solving a system of linear equations</font>
  27. <a name="line18"> 18: </a><font color="#B22222">!T*/</font>
  28. <a name="line19"> 19: </a>! -----------------------------------------------------------------------
  29. <a name="line20"> 20: </a> program main
  30. <a name="line21"> 21: </a> implicit none
  31. <a name="line22"> 22: </a><font color="#A020F0">#include &lt;finclude/petscsys.h&gt;</font>
  32. <a name="line23"> 23: </a><font color="#A020F0">#include &lt;finclude/petscvec.h&gt;</font>
  33. <a name="line24"> 24: </a><font color="#A020F0">#include &lt;finclude/petscmat.h&gt;</font>
  34. <a name="line25"> 25: </a><font color="#A020F0">#include &lt;finclude/petscksp.h&gt;</font>
  35. <a name="line26"> 26: </a><font color="#A020F0">#include &lt;finclude/petscpc.h&gt;</font>
  36. <a name="line27"> 27: </a><font color="#A020F0">#include &lt;finclude/petscviewer.h&gt;</font>
  37. <a name="line28"> 28: </a><font color="#A020F0">#include &lt;finclude/petscviewer.h90&gt;</font>
  38. <a name="line29"> 29: </a> <A href="../../../../../docs/manualpages/Vec/Vec.html#Vec">Vec</A> xvec,bvec,uvec
  39. <a name="line30"> 30: </a> <A href="../../../../../docs/manualpages/Mat/Mat.html#Mat">Mat</A> Amat
  40. <a name="line31"> 31: </a> <A href="../../../../../docs/manualpages/KSP/KSP.html#KSP">KSP</A> ksp
  41. <a name="line32"> 32: </a> <A href="../../../../../docs/manualpages/PC/PC.html#PC">PC</A> pc
  42. <a name="line33"> 33: </a> <A href="../../../../../docs/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</A> ierr
  43. <a name="line34"> 34: </a> <A href="../../../../../docs/manualpages/Viewer/PetscViewer.html#PetscViewer">PetscViewer</A> viewer
  44. <a name="line35"> 35: </a> <A href="../../../../../docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</A> qj,qi,ne,M,Istart,Iend,geq,ix
  45. <a name="line36"> 36: </a> <A href="../../../../../docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</A> ki,kj,lint,nel,ll,j1,i1,ndf
  46. <a name="line37"> 37: </a> <A href="../../../../../docs/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</A> :: idx(4)
  47. <a name="line38"> 38: </a> <A href="../../../../../docs/manualpages/Sys/PetscBool.html#PetscBool">PetscBool</A> flg,out_matlab
  48. <a name="line39"> 39: </a> <A href="../../../../../docs/manualpages/Sys/PetscMPIInt.html#PetscMPIInt">PetscMPIInt</A> npe,mype
  49. <a name="line40"> 40: </a><strong><font color="#FF0000"> <A href="../../../../../docs/manualpages/Sys/PetscScalar.html#PetscScalar">PetscScalar</A>:</font></strong>:ss(4,4),res(4),val
  50. <a name="line41"> 41: </a><strong><font color="#FF0000"> <A href="../../../../../docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</A>:</font></strong>:shp(3,9),sg(3,9),flux(2)
  51. <a name="line42"> 42: </a><strong><font color="#FF0000"> <A href="../../../../../docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</A>:</font></strong>:thk,a1,a2
  52. <a name="line43"> 43: </a> <A href="../../../../../docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</A>, external :: ex54_psi
  53. <a name="line44"> 44: </a><strong><font color="#FF0000"> <A href="../../../../../docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</A>:</font></strong>:norm,tol,theta,eps,h,x,y,xsj
  54. <a name="line45"> 45: </a><strong><font color="#FF0000"> <A href="../../../../../docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</A>:</font></strong>:coord(2,4),dd(2,2),ev(3),blb(2)
  55. <a name="line46"> 46: </a>
  56. <a name="line47"> 47: </a> common /ex54_theta/ theta
  57. <a name="line48"> 48: </a>! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  58. <a name="line49"> 49: </a>! Beginning of program
  59. <a name="line50"> 50: </a>! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  60. <a name="line51"> 51: </a> call <A href="../../../../../docs/manualpages/Sys/PetscInitialize.html#PetscInitialize">PetscInitialize</A>(PETSC_NULL_CHARACTER,ierr)
  61. <a name="line52"> 52: </a> call <A href="http://www.mcs.anl.gov/mpi/www/www3/MPI_Comm_size.html#MPI_Comm_size">MPI_Comm_size</A>(<A href="../../../../../docs/manualpages/Sys/PETSC_COMM_WORLD.html#PETSC_COMM_WORLD">PETSC_COMM_WORLD</A>,npe,ierr)
  62. <a name="line53"> 53: </a> call <A href="http://www.mcs.anl.gov/mpi/www/www3/MPI_Comm_rank.html#MPI_Comm_rank">MPI_Comm_rank</A>(<A href="../../../../../docs/manualpages/Sys/PETSC_COMM_WORLD.html#PETSC_COMM_WORLD">PETSC_COMM_WORLD</A>,mype,ierr)
  63. <a name="line54"> 54: </a>! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  64. <a name="line55"> 55: </a>! set parameters
  65. <a name="line56"> 56: </a>! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  66. <a name="line57"> 57: </a> ne = 9
  67. <a name="line58"> 58: </a> call <A href="../../../../../docs/manualpages/Sys/PetscOptionsGetInt.html#PetscOptionsGetInt">PetscOptionsGetInt</A>(PETSC_NULL_CHARACTER,'-ne',ne,flg,ierr)
  68. <a name="line59"> 59: </a> h = 2.d0/ne
  69. <a name="line60"> 60: </a> M = (ne+1)*(ne+1)
  70. <a name="line61"> 61: </a> theta = 90.d0
  71. <a name="line62"> 62: </a>! theta is input in degrees
  72. <a name="line63"> 63: </a> call <A href="../../../../../docs/manualpages/Sys/PetscOptionsGetReal.html#PetscOptionsGetReal">PetscOptionsGetReal</A>(PETSC_NULL_CHARACTER,'-theta',theta,flg,
  73. <a name="line64"> 64: </a> $ ierr)
  74. <a name="line65"> 65: </a> theta = theta / 57.2957795
  75. <a name="line66"> 66: </a> eps = 1.d0
  76. <a name="line67"> 67: </a> call <A href="../../../../../docs/manualpages/Sys/PetscOptionsGetReal.html#PetscOptionsGetReal">PetscOptionsGetReal</A>(PETSC_NULL_CHARACTER,'-epsilon',eps,flg,
  77. <a name="line68"> 68: </a> $ ierr)
  78. <a name="line69"> 69: </a> ki = 2
  79. <a name="line70"> 70: </a> call <A href="../../../../../docs/manualpages/Sys/PetscOptionsGetRealArray.html#PetscOptionsGetRealArray">PetscOptionsGetRealArray</A>(PETSC_NULL_CHARACTER,'-blob_center',
  80. <a name="line71"> 71: </a> $ blb,ki,flg,ierr)
  81. <a name="line72"> 72: </a> <font color="#4169E1">if</font>( .not. flg ) then
  82. <a name="line73"> 73: </a> blb(1) = 0.d0
  83. <a name="line74"> 74: </a> blb(2) = 0.d0
  84. <a name="line75"> 75: </a> <font color="#4169E1">else</font> <font color="#4169E1">if</font>( ki .ne. 2 ) then
  85. <a name="line76"> 76: </a> print *, 'error: ', ki,
  86. <a name="line77"> 77: </a> $ ' arguments read for -blob_center. Needs to be two.'
  87. <a name="line78"> 78: </a> endif
  88. <a name="line79"> 79: </a> call <A href="../../../../../docs/manualpages/Sys/PetscOptionsGetBool.html#PetscOptionsGetBool">PetscOptionsGetBool</A>(PETSC_NULL_CHARACTER,'-out_matlab',
  89. <a name="line80"> 80: </a> $ out_matlab,flg,ierr)
  90. <a name="line82"> 82: </a> ev(1) = 1.d0
  91. <a name="line83"> 83: </a> ev(2) = eps*ev(1)
  92. <a name="line84"> 84: </a>! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  93. <a name="line85"> 85: </a>! Compute the matrix and right-hand-side vector that define
  94. <a name="line86"> 86: </a>! the linear system, Ax = b.
  95. <a name="line87"> 87: </a>! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  96. <a name="line88"> 88: </a>! Create matrix. When using <A href="../../../../../docs/manualpages/Mat/MatCreate.html#MatCreate">MatCreate</A>(), the matrix format can
  97. <a name="line89"> 89: </a>! be specified at runtime.
  98. <a name="line90"> 90: </a> call <A href="../../../../../docs/manualpages/Mat/MatCreate.html#MatCreate">MatCreate</A>(<A href="../../../../../docs/manualpages/Sys/PETSC_COMM_WORLD.html#PETSC_COMM_WORLD">PETSC_COMM_WORLD</A>,Amat,ierr)
  99. <a name="line91"> 91: </a> call <A href="../../../../../docs/manualpages/Mat/MatSetSizes.html#MatSetSizes">MatSetSizes</A>( Amat,<A href="../../../../../docs/manualpages/Sys/PETSC_DECIDE.html#PETSC_DECIDE">PETSC_DECIDE</A>, <A href="../../../../../docs/manualpages/Sys/PETSC_DECIDE.html#PETSC_DECIDE">PETSC_DECIDE</A>, M, M, ierr )
  100. <a name="line92"> 92: </a> <font color="#4169E1">if</font>( npe == 1 ) then
  101. <a name="line93"> 93: </a> call <A href="../../../../../docs/manualpages/Mat/MatSetType.html#MatSetType">MatSetType</A>( Amat, <A href="../../../../../docs/manualpages/Mat/MATAIJ.html#MATAIJ">MATAIJ</A>, ierr )
  102. <a name="line94"> 94: </a> <font color="#4169E1">else</font>
  103. <a name="line95"> 95: </a> call <A href="../../../../../docs/manualpages/Mat/MatSetType.html#MatSetType">MatSetType</A>( Amat, <A href="../../../../../docs/manualpages/Mat/MATMPIAIJ.html#MATMPIAIJ">MATMPIAIJ</A>, ierr )
  104. <a name="line96"> 96: </a> endif
  105. <a name="line97"> 97: </a> call <A href="../../../../../docs/manualpages/Mat/MatMPIAIJSetPreallocation.html#MatMPIAIJSetPreallocation">MatMPIAIJSetPreallocation</A>(Amat,9,PETSC_NULL_INTEGER,6,
  106. <a name="line98"> 98: </a> $ PETSC_NULL_INTEGER, ierr)
  107. <a name="line99"> 99: </a> call <A href="../../../../../docs/manualpages/Mat/MatSetFromOptions.html#MatSetFromOptions">MatSetFromOptions</A>( Amat, ierr )
  108. <a name="line100">100: </a> call <A href="../../../../../docs/manualpages/Mat/MatSetUp.html#MatSetUp">MatSetUp</A>( Amat, ierr )
  109. <a name="line101">101: </a> call <A href="../../../../../docs/manualpages/Mat/MatGetOwnershipRange.html#MatGetOwnershipRange">MatGetOwnershipRange</A>( Amat, Istart, Iend, ierr )
  110. <a name="line102">102: </a>! Create vectors. Note that we form 1 vector from scratch and
  111. <a name="line103">103: </a>! then duplicate as needed.
  112. <a name="line104">104: </a> call <A href="../../../../../docs/manualpages/Mat/MatGetVecs.html#MatGetVecs">MatGetVecs</A>( Amat, PETSC_NULL_OBJECT, xvec, ierr )
  113. <a name="line105">105: </a> call <A href="../../../../../docs/manualpages/Vec/VecSetFromOptions.html#VecSetFromOptions">VecSetFromOptions</A>( xvec, ierr )
  114. <a name="line106">106: </a> call <A href="../../../../../docs/manualpages/Vec/VecDuplicate.html#VecDuplicate">VecDuplicate</A>( xvec, bvec, ierr )
  115. <a name="line107">107: </a> call <A href="../../../../../docs/manualpages/Vec/VecDuplicate.html#VecDuplicate">VecDuplicate</A>( xvec, uvec, ierr )
  116. <a name="line108">108: </a>! Assemble matrix.
  117. <a name="line109">109: </a>! - Note that <A href="../../../../../docs/manualpages/Mat/MatSetValues.html#MatSetValues">MatSetValues</A>() uses 0-based row and column numbers
  118. <a name="line110">110: </a>! in Fortran as well as in C (as set here in the array <font color="#666666">"col"</font>).
  119. <a name="line111">111: </a> thk = 1.d0 ! thickness
  120. <a name="line112">112: </a> nel = 4 ! nodes per element (quad)
  121. <a name="line113">113: </a> ndf = 1
  122. <a name="line114">114: </a> call int2d(2,lint,sg)
  123. <a name="line115">115: </a> ix = 0
  124. <a name="line116">116: </a> <font color="#4169E1">do</font> geq=Istart,Iend-1,1
  125. <a name="line117">117: </a> qj = geq/(ne+1); qi = mod(geq,(ne+1))
  126. <a name="line118">118: </a> x = h*qi - 1.d0; y = h*qj - 1.d0 ! lower left corner (-1,-1)
  127. <a name="line119">119: </a> <font color="#4169E1">if</font>( qi &lt; ne .and. qj &lt; ne ) then
  128. <a name="line120">120: </a> coord(1,1) = x; coord(2,1) = y
  129. <a name="line121">121: </a> coord(1,2) = x+h; coord(2,2) = y
  130. <a name="line122">122: </a> coord(1,3) = x+h; coord(2,3) = y+h
  131. <a name="line123">123: </a> coord(1,4) = x; coord(2,4) = y+h
  132. <a name="line124">124: </a>! form stiff
  133. <a name="line125">125: </a> ss = 0.d0
  134. <a name="line126">126: </a> <font color="#4169E1">do</font> ll = 1,lint
  135. <a name="line127">127: </a> call shp2dquad(sg(1,ll),sg(2,ll),coord,shp,xsj,2)
  136. <a name="line128">128: </a> xsj = xsj*sg(3,ll)*thk
  137. <a name="line129">129: </a> call thfx2d(ev,coord,shp,dd,2,2,4,ex54_psi)
  138. <a name="line130">130: </a> j1 = 1
  139. <a name="line131">131: </a> <font color="#4169E1">do</font> kj = 1,nel
  140. <a name="line132">132: </a> a1 = (dd(1,1)*shp(1,kj) + dd(1,2)*shp(2,kj))*xsj
  141. <a name="line133">133: </a> a2 = (dd(2,1)*shp(1,kj) + dd(2,2)*shp(2,kj))*xsj
  142. <a name="line134">134: </a><strong><font color="#4169E1"><a name="p"></a>c Compute residual </font></strong>
  143. <a name="line135">135: </a><strong><font color="#4169E1">c p(j1)</font></strong> = p(j1) - a1*gradt(1) - a2*gradt(2)
  144. <a name="line136">136: </a>c Compute tangent
  145. <a name="line137">137: </a> i1 = 1
  146. <a name="line138">138: </a> <font color="#4169E1">do</font> ki = 1,nel
  147. <a name="line139">139: </a> ss(i1,j1) = ss(i1,j1) + a1*shp(1,ki) + a2*shp(2,ki)
  148. <a name="line140">140: </a> i1 = i1 + ndf
  149. <a name="line141">141: </a> end <font color="#4169E1">do</font>
  150. <a name="line142">142: </a> j1 = j1 + ndf
  151. <a name="line143">143: </a> end <font color="#4169E1">do</font>
  152. <a name="line144">144: </a> enddo
  153. <a name="line146">146: </a> idx(1) = geq; idx(2) = geq+1; idx(3) = geq+(ne+1)+1
  154. <a name="line147">147: </a> idx(4) = geq+(ne+1)
  155. <a name="line148">148: </a> <font color="#4169E1">if</font>( qj &gt; 0 ) then
  156. <a name="line149">149: </a> call <A href="../../../../../docs/manualpages/Mat/MatSetValues.html#MatSetValues">MatSetValues</A>(Amat,4,idx,4,idx,ss,<A href="../../../../../docs/manualpages/Sys/ADD_VALUES.html#ADD_VALUES">ADD_VALUES</A>,ierr)
  157. <a name="line150">150: </a> <font color="#4169E1">else</font> ! a BC
  158. <a name="line151">151: </a> <font color="#4169E1">do</font> ki=1,4,1
  159. <a name="line152">152: </a> <font color="#4169E1">do</font> kj=1,4,1
  160. <a name="line153">153: </a> <font color="#4169E1">if</font>(ki&lt;3 .or. kj&lt;3 ) then
  161. <a name="line154">154: </a> <font color="#4169E1">if</font>( ki==kj ) then
  162. <a name="line155">155: </a> ss(ki,kj) = .1d0*ss(ki,kj)
  163. <a name="line156">156: </a> <font color="#4169E1">else</font>
  164. <a name="line157">157: </a> ss(ki,kj) = 0.d0
  165. <a name="line158">158: </a> endif
  166. <a name="line159">159: </a> endif
  167. <a name="line160">160: </a> enddo
  168. <a name="line161">161: </a> enddo
  169. <a name="line162">162: </a> call <A href="../../../../../docs/manualpages/Mat/MatSetValues.html#MatSetValues">MatSetValues</A>(Amat,4,idx,4,idx,ss,<A href="../../../../../docs/manualpages/Sys/ADD_VALUES.html#ADD_VALUES">ADD_VALUES</A>,ierr)
  170. <a name="line163">163: </a> endif ! BC
  171. <a name="line164">164: </a> endif ! add element
  172. <a name="line165">165: </a> <font color="#4169E1">if</font>( qj &gt; 0 ) then ! set rhs
  173. <a name="line167">167: </a> val = h*h*exp(-1.d2*((x+h/2)-blb(1))**2)*
  174. <a name="line168">168: </a> $ exp(-1.d2*((y+h/2)-blb(2))**2)
  175. <a name="line169">169: </a> call <A href="../../../../../docs/manualpages/Vec/VecSetValues.html#VecSetValues">VecSetValues</A>(bvec,1,geq,val,<A href="../../../../../docs/manualpages/Sys/INSERT_VALUES.html#INSERT_VALUES">INSERT_VALUES</A>,ierr)
  176. <a name="line170">170: </a> endif
  177. <a name="line171">171: </a> enddo
  178. <a name="line172">172: </a> call <A href="../../../../../docs/manualpages/Mat/MatAssemblyBegin.html#MatAssemblyBegin">MatAssemblyBegin</A>(Amat,MAT_FINAL_ASSEMBLY,ierr)
  179. <a name="line173">173: </a> call <A href="../../../../../docs/manualpages/Mat/MatAssemblyEnd.html#MatAssemblyEnd">MatAssemblyEnd</A>(Amat,MAT_FINAL_ASSEMBLY,ierr)
  180. <a name="line174">174: </a> call <A href="../../../../../docs/manualpages/Vec/VecAssemblyBegin.html#VecAssemblyBegin">VecAssemblyBegin</A>(bvec,ierr)
  181. <a name="line175">175: </a> call <A href="../../../../../docs/manualpages/Vec/VecAssemblyEnd.html#VecAssemblyEnd">VecAssemblyEnd</A>(bvec,ierr)
  182. <a name="line177">177: </a>! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  183. <a name="line178">178: </a>! Create the linear solver and set various options
  184. <a name="line179">179: </a>! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  185. <a name="line181">181: </a>! Create linear solver context
  186. <a name="line183">183: </a> call <A href="../../../../../docs/manualpages/KSP/KSPCreate.html#KSPCreate">KSPCreate</A>(<A href="../../../../../docs/manualpages/Sys/PETSC_COMM_WORLD.html#PETSC_COMM_WORLD">PETSC_COMM_WORLD</A>,ksp,ierr)
  187. <a name="line185">185: </a>! Set operators. Here the matrix that defines the linear system
  188. <a name="line186">186: </a>! also serves as the preconditioning matrix.
  189. <a name="line188">188: </a> call <A href="../../../../../docs/manualpages/KSP/KSPSetOperators.html#KSPSetOperators">KSPSetOperators</A>(ksp,Amat,Amat,DIFFERENT_NONZERO_PATTERN,ierr)
  190. <a name="line190">190: </a>! Set runtime options, e.g.,
  191. <a name="line191">191: </a>! -ksp_type &lt;type&gt; -pc_type &lt;type&gt; -ksp_monitor -ksp_rtol &lt;rtol&gt;
  192. <a name="line192">192: </a>! These options will override those specified above as long as
  193. <a name="line193">193: </a>! <A href="../../../../../docs/manualpages/KSP/KSPSetFromOptions.html#KSPSetFromOptions">KSPSetFromOptions</A>() is called _after_ any other customization
  194. <a name="line194">194: </a>! routines.
  195. <a name="line196">196: </a> call <A href="../../../../../docs/manualpages/KSP/KSPSetFromOptions.html#KSPSetFromOptions">KSPSetFromOptions</A>(ksp,ierr)
  196. <a name="line198">198: </a>! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  197. <a name="line199">199: </a>! Solve the linear system
  198. <a name="line200">200: </a>! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  199. <a name="line202">202: </a> call <A href="../../../../../docs/manualpages/KSP/KSPSolve.html#KSPSolve">KSPSolve</A>(ksp,bvec,xvec,ierr)
  200. <a name="line205">205: </a>! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  201. <a name="line206">206: </a>! output
  202. <a name="line207">207: </a>! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  203. <a name="line208">208: </a> <font color="#4169E1">if</font>( out_matlab ) then
  204. <a name="line209">209: </a> call <A href="../../../../../docs/manualpages/Viewer/PetscViewerBinaryOpen.html#PetscViewerBinaryOpen">PetscViewerBinaryOpen</A>(<A href="../../../../../docs/manualpages/Sys/PETSC_COMM_WORLD.html#PETSC_COMM_WORLD">PETSC_COMM_WORLD</A>,'Amat',
  205. <a name="line210">210: </a> $ FILE_MODE_WRITE,viewer,ierr)
  206. <a name="line211">211: </a> call <A href="../../../../../docs/manualpages/Mat/MatView.html#MatView">MatView</A>(Amat,viewer,ierr)
  207. <a name="line212">212: </a> call <A href="../../../../../docs/manualpages/Viewer/PetscViewerDestroy.html#PetscViewerDestroy">PetscViewerDestroy</A>(viewer,ierr)
  208. <a name="line214">214: </a> call <A href="../../../../../docs/manualpages/Viewer/PetscViewerBinaryOpen.html#PetscViewerBinaryOpen">PetscViewerBinaryOpen</A>(<A href="../../../../../docs/manualpages/Sys/PETSC_COMM_WORLD.html#PETSC_COMM_WORLD">PETSC_COMM_WORLD</A>,'Bvec',
  209. <a name="line215">215: </a> $ FILE_MODE_WRITE,viewer,ierr)
  210. <a name="line216">216: </a> call <A href="../../../../../docs/manualpages/Vec/VecView.html#VecView">VecView</A>(bvec,viewer,ierr)
  211. <a name="line217">217: </a> call <A href="../../../../../docs/manualpages/Viewer/PetscViewerDestroy.html#PetscViewerDestroy">PetscViewerDestroy</A>(viewer,ierr)
  212. <a name="line219">219: </a> call <A href="../../../../../docs/manualpages/Viewer/PetscViewerBinaryOpen.html#PetscViewerBinaryOpen">PetscViewerBinaryOpen</A>(<A href="../../../../../docs/manualpages/Sys/PETSC_COMM_WORLD.html#PETSC_COMM_WORLD">PETSC_COMM_WORLD</A>,'Xvec',
  213. <a name="line220">220: </a> $ FILE_MODE_WRITE,viewer,ierr)
  214. <a name="line221">221: </a> call <A href="../../../../../docs/manualpages/Vec/VecView.html#VecView">VecView</A>(xvec,viewer,ierr)
  215. <a name="line222">222: </a> call <A href="../../../../../docs/manualpages/Viewer/PetscViewerDestroy.html#PetscViewerDestroy">PetscViewerDestroy</A>(viewer,ierr)
  216. <a name="line224">224: </a> call <A href="../../../../../docs/manualpages/Mat/MatMult.html#MatMult">MatMult</A>(Amat,xvec,uvec,ierr)
  217. <a name="line225">225: </a> val = -1.d0
  218. <a name="line226">226: </a> call <A href="../../../../../docs/manualpages/Vec/VecAXPY.html#VecAXPY">VecAXPY</A>(uvec,val,bvec,ierr)
  219. <a name="line227">227: </a> call <A href="../../../../../docs/manualpages/Viewer/PetscViewerBinaryOpen.html#PetscViewerBinaryOpen">PetscViewerBinaryOpen</A>(<A href="../../../../../docs/manualpages/Sys/PETSC_COMM_WORLD.html#PETSC_COMM_WORLD">PETSC_COMM_WORLD</A>,'Rvec',
  220. <a name="line228">228: </a> $ FILE_MODE_WRITE,viewer,ierr)
  221. <a name="line229">229: </a> call <A href="../../../../../docs/manualpages/Vec/VecView.html#VecView">VecView</A>(uvec,viewer,ierr)
  222. <a name="line230">230: </a> call <A href="../../../../../docs/manualpages/Viewer/PetscViewerDestroy.html#PetscViewerDestroy">PetscViewerDestroy</A>(viewer,ierr)
  223. <a name="line232">232: </a> <font color="#4169E1">if</font>( mype == 0 ) then
  224. <a name="line233">233: </a> open(1,file=<font color="#666666">"ex54f.m"</font>, FORM=<font color="#666666">"formatted"</font>)
  225. <a name="line234">234: </a> write (1,*) <font color="#666666">"A = <A href="../../../../../docs/manualpages/Sys/PetscBinaryRead.html#PetscBinaryRead">PetscBinaryRead</A>('Amat');"</font>
  226. <a name="line235">235: </a> write (1,*) <font color="#666666">"[m n] = size(A);"</font>
  227. <a name="line236">236: </a> write (1,*) <font color="#666666">"mm = sqrt(m);"</font>
  228. <a name="line237">237: </a> write (1,*) <font color="#666666">"b = <A href="../../../../../docs/manualpages/Sys/PetscBinaryRead.html#PetscBinaryRead">PetscBinaryRead</A>('Bvec');"</font>
  229. <a name="line238">238: </a> write (1,*) <font color="#666666">"x = <A href="../../../../../docs/manualpages/Sys/PetscBinaryRead.html#PetscBinaryRead">PetscBinaryRead</A>('Xvec');"</font>
  230. <a name="line239">239: </a> write (1,*) <font color="#666666">"r = <A href="../../../../../docs/manualpages/Sys/PetscBinaryRead.html#PetscBinaryRead">PetscBinaryRead</A>('Rvec');"</font>
  231. <a name="line240">240: </a> write (1,*) <font color="#666666">"bb = reshape(b,mm,mm);"</font>
  232. <a name="line241">241: </a> write (1,*) <font color="#666666">"xx = reshape(x,mm,mm);"</font>
  233. <a name="line242">242: </a> write (1,*) <font color="#666666">"rr = reshape(r,mm,mm);"</font>
  234. <a name="line243">243: </a>c write (1,*) <font color="#666666">"imagesc(bb')"</font>
  235. <a name="line244">244: </a>c write (1,*) <font color="#666666">"title('RHS'),"</font>
  236. <a name="line245">245: </a> write (1,*) <font color="#666666">"figure,"</font>
  237. <a name="line246">246: </a> write (1,*) <font color="#666666">"imagesc(xx')"</font>
  238. <a name="line247">247: </a> write (1,2002) eps,theta*57.2957795
  239. <a name="line248">248: </a> write (1,*) <font color="#666666">"figure,"</font>
  240. <a name="line249">249: </a> write (1,*) <font color="#666666">"imagesc(rr')"</font>
  241. <a name="line250">250: </a> write (1,*) <font color="#666666">"title('Residual'),"</font>
  242. <a name="line251">251: </a> close(1)
  243. <a name="line252">252: </a> endif
  244. <a name="line253">253: </a> endif
  245. <a name="line254">254: </a> 2002 format(<font color="#666666">"title('Solution: esp="</font>,d9.3,<font color="#666666">", theta="</font>,g8.3,<font color="#666666">"'),"</font>)
  246. <a name="line255">255: </a>! Free work space. All PETSc objects should be destroyed when they
  247. <a name="line256">256: </a>! are no longer needed.
  248. <a name="line258">258: </a> call <A href="../../../../../docs/manualpages/Vec/VecDestroy.html#VecDestroy">VecDestroy</A>(xvec,ierr)
  249. <a name="line259">259: </a> call <A href="../../../../../docs/manualpages/Vec/VecDestroy.html#VecDestroy">VecDestroy</A>(bvec,ierr)
  250. <a name="line260">260: </a> call <A href="../../../../../docs/manualpages/Vec/VecDestroy.html#VecDestroy">VecDestroy</A>(uvec,ierr)
  251. <a name="line261">261: </a> call <A href="../../../../../docs/manualpages/Mat/MatDestroy.html#MatDestroy">MatDestroy</A>(Amat,ierr)
  252. <a name="line262">262: </a> call <A href="../../../../../docs/manualpages/KSP/KSPDestroy.html#KSPDestroy">KSPDestroy</A>(ksp,ierr)
  253. <a name="line263">263: </a> call <A href="../../../../../docs/manualpages/Sys/PetscFinalize.html#PetscFinalize">PetscFinalize</A>(ierr)
  254. <a name="line265">265: </a> end
  255. <a name="line267">267: </a>! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  256. <a name="line268">268: </a>! thfx2d - compute material tensor
  257. <a name="line269">269: </a>! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  258. <a name="line271">271: </a> subroutine thfx2d(ev,xl,shp,dd,ndm,ndf,nel,dir)
  259. <a name="line272">272: </a>
  260. <a name="line273">273: </a>c Compute thermal gradient and flux
  261. <a name="line274">274: </a>
  262. <a name="line275">275: </a> implicit none
  263. <a name="line277">277: </a> integer ndm,ndf,nel,i
  264. <a name="line278">278: </a> <A href="../../../../../docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</A> ev(2),xl(ndm,nel),shp(3,*),dir
  265. <a name="line279">279: </a> <A href="../../../../../docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</A> xx,yy,psi,cs,sn,c2,s2,dd(2,2)
  266. <a name="line281">281: </a>c temp = 0.0d0
  267. <a name="line282">282: </a>c gradt(1) = 0.0d0
  268. <a name="line283">283: </a>c gradt(2) = 0.0d0
  269. <a name="line284">284: </a> xx = 0.0d0
  270. <a name="line285">285: </a> yy = 0.0d0
  271. <a name="line286">286: </a> <font color="#4169E1">do</font> i = 1,nel
  272. <a name="line287">287: </a>c gradt(1) = gradt(1) + shp(1,i)*ul(1,i)
  273. <a name="line288">288: </a>c gradt(2) = gradt(2) + shp(2,i)*ul(1,i)
  274. <a name="line289">289: </a>c temp = temp + shp(3,i)*ul(1,i)
  275. <a name="line290">290: </a> xx = xx + shp(3,i)*xl(1,i)
  276. <a name="line291">291: </a> yy = yy + shp(3,i)*xl(2,i)
  277. <a name="line292">292: </a> end <font color="#4169E1">do</font>
  278. <a name="line293">293: </a> psi = dir(xx,yy)
  279. <a name="line294">294: </a>c Compute thermal flux
  280. <a name="line295">295: </a> cs = cos(psi)
  281. <a name="line296">296: </a> sn = sin(psi)
  282. <a name="line297">297: </a> c2 = cs*cs
  283. <a name="line298">298: </a> s2 = sn*sn
  284. <a name="line299">299: </a> cs = cs*sn
  285. <a name="line301">301: </a> dd(1,1) = c2*ev(1) + s2*ev(2)
  286. <a name="line302">302: </a> dd(2,2) = s2*ev(1) + c2*ev(2)
  287. <a name="line303">303: </a> dd(1,2) = cs*(ev(1) - ev(2))
  288. <a name="line304">304: </a> dd(2,1) = dd(1,2)
  289. <a name="line306">306: </a>c flux(1) = -dd(1,1)*gradt(1) - dd(1,2)*gradt(2)
  290. <a name="line307">307: </a>c flux(2) = -dd(2,1)*gradt(1) - dd(2,2)*gradt(2)
  291. <a name="line309">309: </a> end
  292. <a name="line311">311: </a>! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  293. <a name="line312">312: </a>! shp2dquad - shape functions - compute derivatives w/r natural coords.
  294. <a name="line313">313: </a>! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  295. <a name="line314">314: </a> subroutine shp2dquad(s,t,xl,shp,xsj,ndm)
  296. <a name="line315">315: </a>c-----[--.----+----.----+----.-----------------------------------------]
  297. <a name="line316">316: </a>c Purpose: Shape function routine <font color="#4169E1">for</font> 4-node isoparametric quads
  298. <a name="line317">317: </a>c
  299. <a name="line318">318: </a>c Inputs:
  300. <a name="line319">319: </a>c s,t - Natural coordinates of point
  301. <a name="line320">320: </a><strong><font color="#4169E1"><a name="xl"></a>c xl(ndm,*)</font></strong> - Nodal coordinates <font color="#4169E1">for</font> element
  302. <a name="line321">321: </a>c ndm - Spatial dimension of mesh
  303. <a name="line323">323: </a>c Outputs:
  304. <a name="line324">324: </a>c shp(3,*) - Shape functions and derivatives at point
  305. <a name="line325">325: </a>c shp(1,i) = dN_i/dx or dN_i/dxi_1
  306. <a name="line326">326: </a>c shp(2,i) = dN_i/dy or dN_i/dxi_2
  307. <a name="line327">327: </a>c shp(3,i) = N_i
  308. <a name="line328">328: </a>c xsj - Jacobian determinant at point
  309. <a name="line329">329: </a>c-----[--.----+----.----+----.-----------------------------------------]
  310. <a name="line330">330: </a> implicit none
  311. <a name="line331">331: </a> integer ndm
  312. <a name="line332">332: </a> real*8 xo,xs,xt, yo,ys,yt, xsm,xsp,xtm,xtp, ysm,ysp,ytm,ytp
  313. <a name="line333">333: </a> real*8 s,t, xsj,xsj1, sh,th,sp,tp,sm,tm, xl(ndm,4),shp(3,4)
  314. <a name="line335">335: </a>c Set up interpolations
  315. <a name="line337">337: </a> sh = 0.5d0*s
  316. <a name="line338">338: </a> th = 0.5d0*t
  317. <a name="line339">339: </a> sp = 0.5d0 + sh
  318. <a name="line340">340: </a> tp = 0.5d0 + th
  319. <a name="line341">341: </a> sm = 0.5d0 - sh
  320. <a name="line342">342: </a> tm = 0.5d0 - th
  321. <a name="line343">343: </a> shp(3,1) = sm*tm
  322. <a name="line344">344: </a> shp(3,2) = sp*tm
  323. <a name="line345">345: </a> shp(3,3) = sp*tp
  324. <a name="line346">346: </a> shp(3,4) = sm*tp
  325. <a name="line348">348: </a>c Set up natural coordinate functions (times 4)
  326. <a name="line350">350: </a> xo = xl(1,1)-xl(1,2)+xl(1,3)-xl(1,4)
  327. <a name="line351">351: </a> xs = -xl(1,1)+xl(1,2)+xl(1,3)-xl(1,4) + xo*t
  328. <a name="line352">352: </a> xt = -xl(1,1)-xl(1,2)+xl(1,3)+xl(1,4) + xo*s
  329. <a name="line353">353: </a> yo = xl(2,1)-xl(2,2)+xl(2,3)-xl(2,4)
  330. <a name="line354">354: </a> ys = -xl(2,1)+xl(2,2)+xl(2,3)-xl(2,4) + yo*t
  331. <a name="line355">355: </a> yt = -xl(2,1)-xl(2,2)+xl(2,3)+xl(2,4) + yo*s
  332. <a name="line357">357: </a>c Compute jacobian (times 16)
  333. <a name="line359">359: </a> xsj1 = xs*yt - xt*ys
  334. <a name="line361">361: </a>c Divide jacobian by 16 (multiply by .0625)
  335. <a name="line363">363: </a> xsj = 0.0625d0*xsj1
  336. <a name="line364">364: </a> <font color="#4169E1">if</font>(xsj1.eq.0.0d0) then
  337. <a name="line365">365: </a> xsj1 = 1.0d0
  338. <a name="line366">366: </a> <font color="#4169E1">else</font>
  339. <a name="line367">367: </a> xsj1 = 1.0d0/xsj1
  340. <a name="line368">368: </a> endif
  341. <a name="line369">369: </a>
  342. <a name="line370">370: </a>c Divide functions by jacobian
  343. <a name="line371">371: </a>
  344. <a name="line372">372: </a> xs = (xs+xs)*xsj1
  345. <a name="line373">373: </a> xt = (xt+xt)*xsj1
  346. <a name="line374">374: </a> ys = (ys+ys)*xsj1
  347. <a name="line375">375: </a> yt = (yt+yt)*xsj1
  348. <a name="line376">376: </a>
  349. <a name="line377">377: </a>c Multiply by interpolations
  350. <a name="line378">378: </a>
  351. <a name="line379">379: </a> ytm = yt*tm
  352. <a name="line380">380: </a> ysm = ys*sm
  353. <a name="line381">381: </a> ytp = yt*tp
  354. <a name="line382">382: </a> ysp = ys*sp
  355. <a name="line383">383: </a> xtm = xt*tm
  356. <a name="line384">384: </a> xsm = xs*sm
  357. <a name="line385">385: </a> xtp = xt*tp
  358. <a name="line386">386: </a> xsp = xs*sp
  359. <a name="line387">387: </a>
  360. <a name="line388">388: </a>c Compute shape functions
  361. <a name="line389">389: </a>
  362. <a name="line390">390: </a> shp(1,1) = - ytm+ysm
  363. <a name="line391">391: </a> shp(1,2) = ytm+ysp
  364. <a name="line392">392: </a> shp(1,3) = ytp-ysp
  365. <a name="line393">393: </a> shp(1,4) = - ytp-ysm
  366. <a name="line394">394: </a> shp(2,1) = xtm-xsm
  367. <a name="line395">395: </a> shp(2,2) = - xtm-xsp
  368. <a name="line396">396: </a> shp(2,3) = - xtp+xsp
  369. <a name="line397">397: </a> shp(2,4) = xtp+xsm
  370. <a name="line398">398: </a>
  371. <a name="line399">399: </a> end
  372. <a name="line401">401: </a>! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  373. <a name="line402">402: </a>! int2d
  374. <a name="line403">403: </a>! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  375. <a name="line404">404: </a> subroutine int2d(l,lint,sg)
  376. <a name="line405">405: </a>c-----[--.----+----.----+----.-----------------------------------------]
  377. <a name="line406">406: </a>c Purpose: Form Gauss points and weights <font color="#4169E1">for</font> two dimensions
  378. <a name="line407">407: </a>
  379. <a name="line408">408: </a>c Inputs:
  380. <a name="line409">409: </a>c l - Number of points/direction
  381. <a name="line410">410: </a>
  382. <a name="line411">411: </a>c Outputs:
  383. <a name="line412">412: </a>c lint - Total number of points
  384. <a name="line413">413: </a>c sg(3,*) - Array of points and weights
  385. <a name="line414">414: </a>c-----[--.----+----.----+----.-----------------------------------------]
  386. <a name="line415">415: </a> implicit none
  387. <a name="line416">416: </a> integer l,i,lint,lr(9),lz(9)
  388. <a name="line417">417: </a> real*8 g,third,sg(3,*)
  389. <a name="line418">418: </a> data lr/-1,1,1,-1,0,1,0,-1,0/,lz/-1,-1,1,1,-1,0,1,0,0/
  390. <a name="line419">419: </a> data third / 0.3333333333333333d0 /
  391. <a name="line420">420: </a>
  392. <a name="line421">421: </a>c Set number of total points
  393. <a name="line422">422: </a>
  394. <a name="line423">423: </a> lint = l*l
  395. <a name="line424">424: </a>
  396. <a name="line425">425: </a>c 2x2 integration
  397. <a name="line426">426: </a> g = sqrt(third)
  398. <a name="line427">427: </a> <font color="#4169E1">do</font> i = 1,4
  399. <a name="line428">428: </a> sg(1,i) = g*lr(i)
  400. <a name="line429">429: </a> sg(2,i) = g*lz(i)
  401. <a name="line430">430: </a> sg(3,i) = 1.d0
  402. <a name="line431">431: </a> end <font color="#4169E1">do</font>
  403. <a name="line433">433: </a> end
  404. <a name="line435">435: </a>! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  405. <a name="line436">436: </a>! ex54_psi - anusotropic material direction
  406. <a name="line437">437: </a>! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  407. <a name="line438">438: </a> <A href="../../../../../docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</A> function ex54_psi(x,y)
  408. <a name="line439">439: </a> implicit none
  409. <a name="line440">440: </a> <A href="../../../../../docs/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</A> x,y,theta
  410. <a name="line441">441: </a> common /ex54_theta/ theta
  411. <a name="line442">442: </a> ex54_psi = theta
  412. <a name="line443">443: </a> <font color="#4169E1">if</font>( theta &lt; 0. ) then ! circular
  413. <a name="line444">444: </a> <font color="#4169E1">if</font>(y==0) then
  414. <a name="line445">445: </a> ex54_psi = 2.d0*atan(1.d0)
  415. <a name="line446">446: </a> <font color="#4169E1">else</font>
  416. <a name="line447">447: </a> ex54_psi = atan(-x/y)
  417. <a name="line448">448: </a> endif
  418. <a name="line449">449: </a> endif
  419. <a name="line450">450: </a> end
  420. </pre>
  421. </body>
  422. </html>