Actual source code: plate2f.F
1: ! "$Id$";
2: !
3: ! Program usage: mpirun -np <proc> plate2f [all TAO options]
4: !
5: ! This example demonstrates use of the TAO package to solve a bound constrained
6: ! minimization problem. This example is based on a problem from the
7: ! MINPACK-2 test suite. Given a rectangular 2-D domain and boundary values
8: ! along the edges of the domain, the objective is to find the surface
9: ! with the minimal area that satisfies the boundary conditions.
10: ! The command line options are:
11: ! -mx <xg>, where <xg> = number of grid points in the 1st coordinate direction
12: ! -my <yg>, where <yg> = number of grid points in the 2nd coordinate direction
13: ! -bmx <bxg>, where <bxg> = number of grid points under plate in 1st direction
14: ! -bmy <byg>, where <byg> = number of grid points under plate in 2nd direction
15: ! -bheight <ht>, where <ht> = height of the plate
16: !
17: !/*T
18: ! Concepts: TAO - Solving a bound constrained minimization problem
19: ! Routines: TaoInitialize(); TaoFinalize();
20: ! Routines: TaoCreate(); TaoDestroy();
21: ! Routines: TaoAppSetObjectiveAndGradientRoutine();
22: ! Routines: TaoAppSetHessianMat(); TaoAppSetHessianRoutine();
23: ! Routines: TaoAppSetInitialSolutionVec(); TaoAppSetVariableBounds();
24: ! Routines: TaoSetOptions();
25: ! Routines: TaoApplicationCreate(); TaoSolve();
26: ! Routines: TaoView(); TaoAppDestroy();
27: ! Processors: n
28: !T*/
32: implicit none
34: #include "plate2f.h"
36: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
37: ! Variable declarations
38: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
39: !
40: ! Variables:
41: ! (common from plate2f.h):
42: ! Nx, Ny number of processors in x- and y- directions
43: ! mx, my number of grid points in x,y directions
44: ! N global dimension of vector
46: integer info ! used to check for functions returning nonzeros
47: Vec x ! solution vector
48: Vec xl, xu ! lower and upper bounds vectorsp
49: integer m ! number of local elements in vector
50: TAO_SOLVER tao ! TAO_SOLVER solver context
51: TAO_APPLICATION plateapp ! PETSc application
52: Mat H ! Hessian matrix
53: ISLocalToGlobalMapping isltog ! local to global mapping object
54: PetscTruth flg
57: external FormFunctionGradient
58: external FormHessian
59: external MSA_BoundaryConditions
60: external MSA_Plate
61: external MSA_InitialPoint
62: ! Initialize Tao
63: call PetscInitialize(PETSC_NULL_CHARACTER,info)
64: call TaoInitialize(PETSC_NULL_CHARACTER,info)
65:
67: ! Specify default dimensions of the problem
68: mx = 10
69: my = 10
70: bheight = 0.1
72: ! Check for any command line arguments that override defaults
73:
74: call PetscOptionsGetInt(TAO_NULL_CHARACTER,"-mx",mx,flg,info)
75: call PetscOptionsGetInt(TAO_NULL_CHARACTER,"-my",my,flg,info)
76:
77: bmx = mx/2
78: bmy = my/2
80: call PetscOptionsGetInt(TAO_NULL_CHARACTER,"-bmx",bmx,flg,info)
81: call PetscOptionsGetInt(TAO_NULL_CHARACTER,"-bmy",bmy,flg,info)
82: call PetscOptionsGetReal(TAO_NULL_CHARACTER,"-bheight",bheight, &
83: & flg,info)
84:
86: ! Calculate any derived values from parameters
87: N = mx*my
89: ! Let Petsc determine the dimensions of the local vectors
90: Nx = PETSC_DECIDE
91: NY = PETSC_DECIDE
93: ! A two dimensional distributed array will help define this problem, which
94: ! derives from an elliptic PDE on a two-dimensional domain. From the
95: ! distributed array, create the vectors
97: call DACreate2d(MPI_COMM_WORLD,DA_NONPERIODIC,DA_STENCIL_BOX, &
98: & mx,my,Nx,Ny,1,1,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER, &
99: & da,info)
100:
102: ! Extract global and local vectors from DA; The local vectors are
103: ! used solely as work space for the evaluation of the function,
104: ! gradient, and Hessian. Duplicate for remaining vectors that are
105: ! the same types.
107: call DACreateGlobalVector(da,x,info)
108: call DACreateLocalVector(da,localX,info)
109: call VecDuplicate(localX,localV,info)
111: ! Create a matrix data structure to store the Hessian.
112: ! Here we (optionally) also associate the local numbering scheme
113: ! with the matrix so that later we can use local indices for matrix
114: ! assembly
116: call VecGetLocalSize(x,m,info)
117: call MatCreateMPIAIJ(MPI_COMM_WORLD,m,m,N,N,7,PETSC_NULL_INTEGER, &
118: & 3,PETSC_NULL_INTEGER,H,info)
120: call MatSetOption(H,MAT_SYMMETRIC,info)
121: call DAGetISLocalToGlobalMapping(da,isltog,info)
122: call MatSetLocalToGlobalMapping(H,isltog,info)
123:
125: ! The Tao code begins here
126: ! Create TAO solver and set desired solution method.
127: ! This problems uses bounded variables, so the
128: ! method must either be 'tao_tron' or 'tao_blmvm'
130: call TaoCreate(MPI_COMM_WORLD,'tao_blmvm',tao,info)
131: call TaoApplicationCreate(MPI_COMM_WORLD,plateapp,info)
132:
134: ! Set minimization function and gradient, hessian evaluation functions
136: ! TaoAppSetObjectiveAndGradientRoutine is shortened to 31 chars to comply with some compilers
137: call TaoAppSetObjectiveAndGradientRo(plateapp, &
138: & FormFunctionGradient,PETSC_NULL_OBJECT,info)
140: call TaoAppSetHessianMat(plateapp,H,H,info)
141: call TaoAppSetHessianRoutine(plateapp,FormHessian, &
142: & PETSC_NULL_OBJECT, info)
143:
145: ! Set Variable bounds
146: call MSA_BoundaryConditions(info)
147: call VecDuplicate(x,xl,info)
148: call VecDuplicate(x,xu,info)
149: call MSA_Plate(xl,xu,info)
150: call TaoAppSetVariableBounds(plateapp,xl,xu,info)
152: ! Set the initial solution guess
153: call MSA_InitialPoint(x, info)
154: call TaoAppSetInitialSolutionVec(plateapp,x,info)
156: ! Check for any tao command line options
157: call TaoSetOptions(plateapp,tao,info)
159: ! Solve the application
160: call TaoSolveApplication(plateapp,tao,info)
162: ! View TAO solver information
163: ! call TaoView(tao,info)
165: ! Free TAO data structures
166: call TaoDestroy(tao,info)
167: call TaoAppDestroy(plateapp,info)
169: ! Free PETSc data structures
170: call VecDestroy(x,info)
171: call VecDestroy(xl,info)
172: call VecDestroy(xu,info)
173: call VecDestroy(Top,info)
174: call VecDestroy(Bottom,info)
175: call VecDestroy(Left,info)
176: call VecDestroy(Right,info)
177: call MatDestroy(H,info)
178: call VecDestroy(localX,info)
179: call VecDestroy(localV,info)
180: call DADestroy(da,info)
182: ! Finalize TAO
184: call TaoFinalize(info)
185: call PetscFinalize(info)
187: end
189: ! ---------------------------------------------------------------------
190: !
191: ! FormFunctionGradient - Evaluates function f(X).
192: !
193: ! Input Parameters:
194: ! tao - the TAO_SOLVER context
195: ! X - the input vector
196: ! dummy - optional user-defined context, as set by TaoSetFunction()
197: ! (not used here)
198: !
199: ! Output Parameters:
200: ! fcn - the newly evaluated function
201: ! G - the gradient vector
202: ! info - error code
203: !
206: subroutine FormFunctionGradient(tao,X,fcn,G,dummy,info)
207: implicit none
209: ! da, localX, localG, Top, Bottom, Left, Right defined in plate2f.h
210: #include "plate2f.h"
211:
212: ! Input/output variables
214: TAO_SOLVER tao
215: PetscScalar fcn
216: Vec X, G
217: integer dummy, info
218:
219: integer i,j,row
220: integer xs, xm, gxs, gxm, ys, ym, gys, gym
221: PetscScalar ft,zero,hx,hy,hydhx,hxdhy
222: PetscScalar area,rhx,rhy
223: PetscScalar f1,f2,f3,f4,f5,f6,d1,d2,d3
224: PetscScalar d4,d5,d6,d7,d8
225: PetscScalar df1dxc,df2dxc,df3dxc,df4dxc
226: PetscScalar df5dxc,df6dxc
227: PetscScalar xc,xl,xr,xt,xb,xlt,xrb
230: ! PETSc's VecGetArray acts differently in Fortran than it does in C.
231: ! Calling VecGetArray((Vec) X, (PetscScalar) x_array(0:1), (PetscOffset) x_index, info)
232: ! will return an array of doubles referenced by x_array offset by x_index.
233: ! i.e., to reference the kth element of X, use x_array(k + x_index).
234: ! Notice that by declaring the arrays with range (0:1), we are using the C 0-indexing practice.
235: PetscScalar g_v(0:1),x_v(0:1)
236: PetscScalar top_v(0:1),left_v(0:1)
237: PetscScalar right_v(0:1),bottom_v(0:1)
238: PetscOffset g_i,left_i,right_i
239: PetscOffset bottom_i,top_i,x_i
241: ft = 0.0d0
242: zero = 0.0d0
243: hx = 1.0d0/(mx + 1)
244: hy = 1.0d0/(my + 1)
245: hydhx = hy/hx
246: hxdhy = hx/hy
247: area = 0.5d0 * hx * hy
248: rhx = mx + 1.0d0
249: rhy = my + 1.0d0
252: ! Get local mesh boundaries
253: call DAGetCorners(da,xs,ys,PETSC_NULL_INTEGER,xm,ym, &
254: & PETSC_NULL_INTEGER,info)
255: call DAGetGhostCorners(da,gxs,gys,PETSC_NULL_INTEGER, &
256: & gxm,gym,PETSC_NULL_INTEGER,info)
258: ! Scatter ghost points to local vector
259: call DAGlobalToLocalBegin(da,X,INSERT_VALUES,localX,info)
260: call DAGlobalToLocalEnd(da,X,INSERT_VALUES,localX,info)
262: ! Initialize the vector to zero
263: call VecSet(localV,zero,info)
265: ! Get arrays to vector data (See note above about using VecGetArray in Fortran)
266: call VecGetArray(localX,x_v,x_i,info)
267: call VecGetArray(localV,g_v,g_i,info)
268: call VecGetArray(Top,top_v,top_i,info)
269: call VecGetArray(Bottom,bottom_v,bottom_i,info)
270: call VecGetArray(Left,left_v,left_i,info)
271: call VecGetArray(Right,right_v,right_i,info)
273: ! Compute function over the locally owned part of the mesh
274: do j = ys,ys+ym-1
275: do i = xs,xs+xm-1
276: row = (j-gys)*gxm + (i-gxs)
277: xc = x_v(row+x_i)
278: xt = xc
279: xb = xc
280: xr = xc
281: xl = xc
282: xrb = xc
283: xlt = xc
285: if (i .eq. 0) then !left side
286: xl = left_v(j - ys + 1 + left_i)
287: xlt = left_v(j - ys + 2 + left_i)
288: else
289: xl = x_v(row - 1 + x_i)
290: endif
292: if (j .eq. 0) then !bottom side
293: xb = bottom_v(i - xs + 1 + bottom_i)
294: xrb = bottom_v(i - xs + 2 + bottom_i)
295: else
296: xb = x_v(row - gxm + x_i)
297: endif
299: if (i + 1 .eq. gxs + gxm) then !right side
300: xr = right_v(j - ys + 1 + right_i)
301: xrb = right_v(j - ys + right_i)
302: else
303: xr = x_v(row + 1 + x_i)
304: endif
306: if (j + 1 .eq. gys + gym) then !top side
307: xt = top_v(i - xs + 1 + top_i)
308: xlt = top_v(i - xs + top_i)
309: else
310: xt = x_v(row + gxm + x_i)
311: endif
313: if ((i .gt. gxs ) .and. (j + 1 .lt. gys + gym)) then
314: xlt = x_v(row - 1 + gxm + x_i)
315: endif
317: if ((j .gt. gys) .and. (i + 1 .lt. gxs + gxm)) then
318: xrb = x_v(row + 1 - gxm + x_i)
319: endif
321: d1 = xc-xl
322: d2 = xc-xr
323: d3 = xc-xt
324: d4 = xc-xb
325: d5 = xr-xrb
326: d6 = xrb-xb
327: d7 = xlt-xl
328: d8 = xt-xlt
330: df1dxc = d1 * hydhx
331: df2dxc = d1 * hydhx + d4 * hxdhy
332: df3dxc = d3 * hxdhy
333: df4dxc = d2 * hydhx + d3 * hxdhy
334: df5dxc = d2 * hydhx
335: df6dxc = d4 * hxdhy
337: d1 = d1 * rhx
338: d2 = d2 * rhx
339: d3 = d3 * rhy
340: d4 = d4 * rhy
341: d5 = d5 * rhy
342: d6 = d6 * rhx
343: d7 = d7 * rhy
344: d8 = d8 * rhx
346: f1 = sqrt(1.0d0 + d1*d1 + d7*d7)
347: f2 = sqrt(1.0d0 + d1*d1 + d4*d4)
348: f3 = sqrt(1.0d0 + d3*d3 + d8*d8)
349: f4 = sqrt(1.0d0 + d3*d3 + d2*d2)
350: f5 = sqrt(1.0d0 + d2*d2 + d5*d5)
351: f6 = sqrt(1.0d0 + d4*d4 + d6*d6)
353: ft = ft + f2 + f4
355: df1dxc = df1dxc / f1
356: df2dxc = df2dxc / f2
357: df3dxc = df3dxc / f3
358: df4dxc = df4dxc / f4
359: df5dxc = df5dxc / f5
360: df6dxc = df6dxc / f6
362: g_v(row + g_i) = 0.5 * (df1dxc + df2dxc + df3dxc + df4dxc + &
363: & df5dxc + df6dxc)
364: enddo
365: enddo
367: ! Compute triangular areas along the border of the domain.
368: if (xs .eq. 0) then ! left side
369: do j=ys,ys+ym-1
370: d3 = (left_v(j-ys+1+left_i) - left_v(j-ys+2+left_i)) &
371: & * rhy
372: d2 = (left_v(j-ys+1+left_i) - x_v((j-gys)*gxm + x_i)) &
373: & * rhx
374: ft = ft + sqrt(1.0d0 + d3*d3 + d2*d2)
375: enddo
376: endif
378:
379: if (ys .eq. 0) then !bottom side
380: do i=xs,xs+xm-1
381: d2 = (bottom_v(i+1-xs+bottom_i)-bottom_v(i-xs+2+bottom_i)) &
382: & * rhx
383: d3 = (bottom_v(i-xs+1+bottom_i)-x_v(i-gxs+x_i))*rhy
384: ft = ft + sqrt(1.0 + d3*d3 + d2*d2)
385: enddo
386: endif
388:
389: if (xs + xm .eq. mx) then ! right side
390: do j=ys,ys+ym-1
391: d1 = (x_v((j+1-gys)*gxm-1+x_i)-right_v(j-ys+1+right_i))*rhx
392: d4 = (right_v(j-ys+right_i) - right_v(j-ys+1+right_i))*rhy
393: ft = ft + sqrt(1.0d0 + d1*d1 + d4*d4)
394: enddo
395: endif
397:
398: if (ys + ym .eq. my) then
399: do i=xs,xs+xm-1
400: d1 = (x_v((gym-1)*gxm+i-gxs+x_i) - top_v(i-xs+1+top_i))*rhy
401: d4 = (top_v(i-xs+1+top_i) - top_v(i-xs+top_i))*rhx
402: ft = ft + sqrt(1.0d0 + d1*d1 + d4*d4)
403: enddo
404: endif
406:
407: if ((ys .eq. 0) .and. (xs .eq. 0)) then
408: d1 = (left_v(0 + left_i) - left_v(1 + left_i)) * rhy
409: d2 = (bottom_v(0+bottom_i)-bottom_v(1+bottom_i))*rhx
410: ft = ft + sqrt(1.0d0 + d1*d1 + d2*d2)
411: endif
413: if ((ys + ym .eq. my) .and. (xs + xm .eq. mx)) then
414: d1 = (right_v(ym+1+right_i) - right_v(ym+right_i))*rhy
415: d2 = (top_v(xm+1+top_i) - top_v(xm + top_i))*rhx
416: ft = ft + sqrt(1.0d0 + d1*d1 + d2*d2)
417: endif
419: ft = ft * area
420: call MPI_Allreduce(ft,fcn,1,MPI_DOUBLE_PRECISION, &
421: & MPI_SUM,MPI_COMM_WORLD,info)
425: ! Restore vectors
426: call VecRestoreArray(localX,x_v,x_i,info)
427: call VecRestoreArray(localV,g_v,g_i,info)
428: call VecRestoreArray(Left,left_v,left_i,info)
429: call VecRestoreArray(Top,top_v,top_i,info)
430: call VecRestoreArray(Bottom,bottom_v,bottom_i,info)
431: call VecRestoreArray(Right,right_v,right_i,info)
433: ! Scatter values to global vector
434: call DALocalToGlobal(da,localV,INSERT_VALUES,G,info)
436: call PetscLogFlops(70*xm*ym,info)
438: return
439: end !FormFunctionGradient
440:
445: ! ----------------------------------------------------------------------------
446: !
447: !/*
448: ! FormHessian - Evaluates Hessian matrix.
449: !
450: ! Input Parameters:
451: !. tao - the TAO_SOLVER context
452: !. X - input vector
453: !. dummy - not used
454: !
455: ! Output Parameters:
456: !. Hessian - Hessian matrix
457: !. Hpc - optionally different preconditioning matrix
458: !. flag - flag indicating matrix structure
459: !
460: ! Notes:
461: ! Due to mesh point reordering with DAs, we must always work
462: ! with the local mesh points, and then transform them to the new
463: ! global numbering with the local-to-global mapping. We cannot work
464: ! directly with the global numbers for the original uniprocessor mesh!
465: !
466: ! Two methods are available for imposing this transformation
467: ! when setting matrix entries:
468: ! (A) MatSetValuesLocal(), using the local ordering (including
469: ! ghost points!)
470: ! - Do the following two steps once, before calling TaoSolve()
471: ! - Use DAGetISLocalToGlobalMapping() to extract the
472: ! local-to-global map from the DA
473: ! - Associate this map with the matrix by calling
474: ! MatSetLocalToGlobalMapping()
475: ! - Then set matrix entries using the local ordering
476: ! by calling MatSetValuesLocal()
477: ! (B) MatSetValues(), using the global ordering
478: ! - Use DAGetGlobalIndices() to extract the local-to-global map
479: ! - Then apply this map explicitly yourself
480: ! - Set matrix entries using the global ordering by calling
481: ! MatSetValues()
482: ! Option (A) seems cleaner/easier in many cases, and is the procedure
483: ! used in this example.
484: */
485: subroutine FormHessian(tao, X, Hessian, Hpc, flg, dummy, info)
486: implicit none
488: ! da,Top,Left,Right,Bottom,mx,my,localX defined in plate2f.h
489: #include "plate2f.h"
490:
491: TAO_SOLVER tao
492: Vec X
493: Mat Hessian,Hpc
494: MatStructure flg
495: integer dummy,info
497: integer i,j,k,row
498: integer xs,xm,gxs,gxm,ys,ym,gys,gym,col(0:6)
499: PetscScalar hx,hy,hydhx,hxdhy,rhx,rhy
500: PetscScalar f1,f2,f3,f4,f5,f6,d1,d2,d3
501: PetscScalar d4,d5,d6,d7,d8
502: PetscScalar xc,xl,xr,xt,xb,xlt,xrb
503: PetscScalar hl,hr,ht,hb,hc,htl,hbr
505: ! PETSc's VecGetArray acts differently in Fortran than it does in C.
506: ! Calling VecGetArray((Vec) X, (PetscScalar) x_array(0:1), (PetscOffset) x_index, info)
507: ! will return an array of doubles referenced by x_array offset by x_index.
508: ! i.e., to reference the kth element of X, use x_array(k + x_index).
509: ! Notice that by declaring the arrays with range (0:1), we are using the C 0-indexing practice.
510: PetscScalar right_v(0:1),left_v(0:1)
511: PetscScalar bottom_v(0:1),top_v(0:1)
512: PetscScalar x_v(0:1)
513: PetscOffset x_i,right_i,left_i
514: PetscOffset bottom_i,top_i
515: PetscScalar v(0:6)
516: PetscTruth assembled
517:
518: ! Set various matrix options
519: call MatSetOption(Hessian,MAT_IGNORE_OFF_PROC_ENTRIES,info)
520: call MatSetOption(Hessian,MAT_COLUMNS_SORTED,info)
521: call MatSetOption(Hessian,MAT_ROWS_SORTED,info)
523: ! Get local mesh boundaries
524: call DAGetCorners(da,xs,ys,PETSC_NULL_INTEGER,xm,ym, &
525: & PETSC_NULL_INTEGER,info)
526: call DAGetGhostCorners(da,gxs,gys,PETSC_NULL_INTEGER,gxm,gym, &
527: & PETSC_NULL_INTEGER,info)
529: ! Scatter ghost points to local vectors
530: call DAGlobalToLocalBegin(da,X,INSERT_VALUES,localX,info)
531: call DAGlobalToLocalEnd(da,X,INSERT_VALUES,localX,info)
533: ! Get pointers to vector data (see note on Fortran arrays above)
534: call VecGetArray(localX,x_v,x_i,info)
535: call VecGetArray(Top,top_v,top_i,info)
536: call VecGetArray(Bottom,bottom_v,bottom_i,info)
537: call VecGetArray(Left,left_v,left_i,info)
538: call VecGetArray(Right,right_v,right_i,info)
540: ! Initialize matrix entries to zero
541: call MatAssembled(Hessian,assembled,info)
542: if (assembled .eq. PETSC_TRUE) call MatZeroEntries(Hessian,info)
545: rhx = mx + 1.0
546: rhy = my + 1.0
547: hx = 1.0/rhx
548: hy = 1.0/rhy
549: hydhx = hy/hx
550: hxdhy = hx/hy
551: ! compute Hessian over the locally owned part of the mesh
553: do i=xs,xs+xm-1
554: do j=ys,ys+ym-1
555: row = (j-gys)*gxm + (i-gxs)
556:
557: xc = x_v(row + x_i)
558: xt = xc
559: xb = xc
560: xr = xc
561: xl = xc
562: xrb = xc
563: xlt = xc
565: if (i .eq. gxs) then ! Left side
566: xl = left_v(left_i + j - ys + 1)
567: xlt = left_v(left_i + j - ys + 2)
568: else
569: xl = x_v(x_i + row -1 )
570: endif
572: if (j .eq. gys) then ! bottom side
573: xb = bottom_v(bottom_i + i - xs + 1)
574: xrb = bottom_v(bottom_i + i - xs + 2)
575: else
576: xb = x_v(x_i + row - gxm)
577: endif
579: if (i+1 .eq. gxs + gxm) then !right side
580: xr = right_v(right_i + j - ys + 1)
581: xrb = right_v(right_i + j - ys)
582: else
583: xr = x_v(x_i + row + 1)
584: endif
586: if (j+1 .eq. gym+gys) then !top side
587: xt = top_v(top_i +i - xs + 1)
588: xlt = top_v(top_i + i - xs)
589: else
590: xt = x_v(x_i + row + gxm)
591: endif
593: if ((i .gt. gxs) .and. (j+1 .lt. gys+gym)) then
594: xlt = x_v(x_i + row - 1 + gxm)
595: endif
596:
597: if ((i+1 .lt. gxs+gxm) .and. (j .gt. gys)) then
598: xrb = x_v(x_i + row + 1 - gxm)
599: endif
601: d1 = (xc-xl)*rhx
602: d2 = (xc-xr)*rhx
603: d3 = (xc-xt)*rhy
604: d4 = (xc-xb)*rhy
605: d5 = (xrb-xr)*rhy
606: d6 = (xrb-xb)*rhx
607: d7 = (xlt-xl)*rhy
608: d8 = (xlt-xt)*rhx
609:
610: f1 = sqrt( 1.0d0 + d1*d1 + d7*d7)
611: f2 = sqrt( 1.0d0 + d1*d1 + d4*d4)
612: f3 = sqrt( 1.0d0 + d3*d3 + d8*d8)
613: f4 = sqrt( 1.0d0 + d3*d3 + d2*d2)
614: f5 = sqrt( 1.0d0 + d2*d2 + d5*d5)
615: f6 = sqrt( 1.0d0 + d4*d4 + d6*d6)
616:
617:
618: hl = (-hydhx*(1.0+d7*d7)+d1*d7)/(f1*f1*f1)+ &
619: & (-hydhx*(1.0+d4*d4)+d1*d4)/(f2*f2*f2)
621: hr = (-hydhx*(1.0+d5*d5)+d2*d5)/(f5*f5*f5)+ &
622: & (-hydhx*(1.0+d3*d3)+d2*d3)/(f4*f4*f4)
624: ht = (-hxdhy*(1.0+d8*d8)+d3*d8)/(f3*f3*f3)+ &
625: & (-hxdhy*(1.0+d2*d2)+d2*d3)/(f4*f4*f4)
627: hb = (-hxdhy*(1.0+d6*d6)+d4*d6)/(f6*f6*f6)+ &
628: & (-hxdhy*(1.0+d1*d1)+d1*d4)/(f2*f2*f2)
629:
630: hbr = -d2*d5/(f5*f5*f5) - d4*d6/(f6*f6*f6)
631: htl = -d1*d7/(f1*f1*f1) - d3*d8/(f3*f3*f3)
632:
633: hc = hydhx*(1.0+d7*d7)/(f1*f1*f1) + &
634: & hxdhy*(1.0+d8*d8)/(f3*f3*f3) + &
635: & hydhx*(1.0+d5*d5)/(f5*f5*f5) + &
636: & hxdhy*(1.0+d6*d6)/(f6*f6*f6) + &
637: & (hxdhy*(1.0+d1*d1)+hydhx*(1.0+d4*d4)- &
638: & 2*d1*d4)/(f2*f2*f2) + (hxdhy*(1.0+d2*d2)+ &
639: & hydhx*(1.0+d3*d3)-2*d2*d3)/(f4*f4*f4)
640:
641: hl = hl * 0.5
642: hr = hr * 0.5
643: ht = ht * 0.5
644: hb = hb * 0.5
645: hbr = hbr * 0.5
646: htl = htl * 0.5
647: hc = hc * 0.5
649: k = 0
651: if (j .gt. 0) then
652: v(k) = hb
653: col(k) = row - gxm
654: k=k+1
655: endif
657: if ((j .gt. 0) .and. (i .lt. mx-1)) then
658: v(k) = hbr
659: col(k) = row-gxm+1
660: k=k+1
661: endif
663: if (i .gt. 0) then
664: v(k) = hl
665: col(k) = row - 1
666: k = k+1
667: endif
669: v(k) = hc
670: col(k) = row
671: k=k+1
673: if (i .lt. mx-1) then
674: v(k) = hr
675: col(k) = row + 1
676: k=k+1
677: endif
679: if ((i .gt. 0) .and. (j .lt. my-1)) then
680: v(k) = htl
681: col(k) = row + gxm - 1
682: k=k+1
683: endif
685: if (j .lt. my-1) then
686: v(k) = ht
687: col(k) = row + gxm
688: k=k+1
689: endif
691: ! Set matrix values using local numbering, defined earlier in main routine
692: call MatSetValuesLocal(Hessian,1,row,k,col,v,INSERT_VALUES, &
693: & info)
695:
697: enddo
698: enddo
699:
700: ! restore vectors
701: call VecRestoreArray(localX,x_v,x_i,info)
702: call VecRestoreArray(Left,left_v,left_i,info)
703: call VecRestoreArray(Right,right_v,right_i,info)
704: call VecRestoreArray(Top,top_v,top_i,info)
705: call VecRestoreArray(Bottom,bottom_v,bottom_i,info)
708: ! Assemble the matrix
709: call MatAssemblyBegin(Hessian,MAT_FINAL_ASSEMBLY,info)
710: call MatAssemblyEnd(Hessian,MAT_FINAL_ASSEMBLY,info)
712: call PetscLogFlops(199*xm*ym,info)
714: return
715: end
716:
717:
721: ! Top,Left,Right,Bottom,bheight,mx,my,bmx,bmy,H, defined in plate2f.h
723: ! ----------------------------------------------------------------------------
724: !
725: !/*
726: ! MSA_BoundaryConditions - calculates the boundary conditions for the region
727: !
728: !
729: !*/
731: subroutine MSA_BoundaryConditions(info)
732: implicit none
734: ! Top,Left,Right,Bottom,bheight,mx,my,bmx,bmy defined in plate2f.h
735: #include "plate2f.h"
737: integer i,j,k,limit,info,maxits
738: integer xs, xm, gxs, gxm, ys, ym, gys, gym
739: integer bsize, lsize, tsize, rsize
740: PetscScalar one,two,three,tol
741: PetscScalar scl,fnorm,det,xt
742: PetscScalar yt,hx,hy,u1,u2,nf1,nf2
743: PetscScalar njac11,njac12,njac21,njac22
744: PetscScalar b, t, l, r
745: PetscScalar boundary_v(0:1)
746: PetscOffset boundary_i
747: logical exitloop
748: TaoTruth flg
750: limit=0
751: maxits = 5
752: tol=1e-10
753: b=-0.5d0
754: t= 0.5d0
755: l=-0.5d0
756: r= 0.5d0
757: xt=0
758: yt=0
759: one=1.0d0
760: two=2.0d0
761: three=3.0d0
764: call DAGetCorners(da,xs,ys,PETSC_NULL_INTEGER,xm,ym, &
765: & PETSC_NULL_INTEGER,info)
766: call DAGetGhostCorners(da,gxs,gys,PETSC_NULL_INTEGER, &
767: & gxm,gym,PETSC_NULL_INTEGER,info)
769: bsize = xm + 2
770: lsize = ym + 2
771: rsize = ym + 2
772: tsize = xm + 2
773:
775: call VecCreateMPI(MPI_COMM_WORLD,bsize,PETSC_DECIDE,Bottom,info)
776: call VecCreateMPI(MPI_COMM_WORLD,tsize,PETSC_DECIDE,Top,info)
777: call VecCreateMPI(MPI_COMM_WORLD,lsize,PETSC_DECIDE,Left,info)
778: call VecCreateMPI(MPI_COMM_WORLD,rsize,PETSC_DECIDE,Right,info)
780: hx= (r-l)/(mx+1)
781: hy= (t-b)/(my+1)
783: do j=0,3
784:
785: if (j.eq.0) then
786: yt=b
787: xt=l+hx*xs
788: limit=bsize
789: call VecGetArray(Bottom,boundary_v,boundary_i,info)
790:
792: elseif (j.eq.1) then
793: yt=t
794: xt=l+hx*xs
795: limit=tsize
796: call VecGetArray(Top,boundary_v,boundary_i,info)
798: elseif (j.eq.2) then
799: yt=b+hy*ys
800: xt=l
801: limit=lsize
802: call VecGetArray(Left,boundary_v,boundary_i,info)
804: elseif (j.eq.3) then
805: yt=b+hy*ys
806: xt=r
807: limit=rsize
808: call VecGetArray(Right,boundary_v,boundary_i,info)
809: endif
810:
812: do i=0,limit-1
813:
814: u1=xt
815: u2=-yt
816: k = 0
817: exitloop = .false.
818: do while (k .lt. maxits .and. (.not. exitloop) )
820: nf1=u1 + u1*u2*u2 - u1*u1*u1/three-xt
821: nf2=-u2 - u1*u1*u2 + u2*u2*u2/three-yt
822: fnorm=sqrt(nf1*nf1+nf2*nf2)
823: if (fnorm .gt. tol) then
824: njac11=one+u2*u2-u1*u1
825: njac12=two*u1*u2
826: njac21=-two*u1*u2
827: njac22=-one - u1*u1 + u2*u2
828: det = njac11*njac22-njac21*njac12
829: u1 = u1-(njac22*nf1-njac12*nf2)/det
830: u2 = u2-(njac11*nf2-njac21*nf1)/det
831: else
832: exitloop = .true.
833: endif
834: k=k+1
835: enddo
837: boundary_v(i + boundary_i) = u1*u1-u2*u2
838: if ((j .eq. 0) .or. (j .eq. 1)) then
839: xt = xt + hx
840: else
841: yt = yt + hy
842: endif
844: enddo
845:
847: if (j.eq.0) then
848: call VecRestoreArray(Bottom,boundary_v,boundary_i,info)
849: elseif (j.eq.1) then
850: call VecRestoreArray(Top,boundary_v,boundary_i,info)
851: elseif (j.eq.2) then
852: call VecRestoreArray(Left,boundary_v,boundary_i,info)
853: elseif (j.eq.3) then
854: call VecRestoreArray(Right,boundary_v,boundary_i,info)
855: endif
856:
857: enddo
860: ! Scale the boundary if desired
861: call PetscOptionsGetReal(PETSC_NULL_CHARACTER,"-bottom", &
862: & scl,flg,info)
863: if (flg .eq. PETSC_TRUE) then
864: call VecScale(scl,Bottom,info)
865: endif
867: call PetscOptionsGetReal(PETSC_NULL_CHARACTER,"-top", &
868: & scl,flg,info)
869: if (flg .eq. PETSC_TRUE) then
870: call VecScale(scl,Top,info)
871: endif
873: call PetscOptionsGetReal(PETSC_NULL_CHARACTER,"-right", &
874: & scl,flg,info)
875: if (flg .eq. PETSC_TRUE) then
876: call VecScale(scl,Right,info)
877: endif
879: call PetscOptionsGetReal(PETSC_NULL_CHARACTER,"-left", &
880: & scl,flg,info)
881: if (flg .eq. PETSC_TRUE) then
882: call VecScale(scl,Left,info)
883: endif
884:
885:
886: return
887: end
889: ! ----------------------------------------------------------------------------
890: !
891: !/*
892: ! MSA_Plate - Calculates an obstacle for surface to stretch over
893: !
894: ! Output Parameter:
895: !. xl - lower bound vector
896: !. xu - upper bound vector
897: !
898: !*/
900: subroutine MSA_Plate(xl,xu,info)
901: implicit none
903: ! mx,my,bmx,bmy,da,bheight defined in plate2f.h
904: #include "plate2f.h"
905: Vec xl,xu
906: integer i,j,row,info
907: integer xs, xm, ys, ym
908: PetscScalar lb,ub
910: ! PETSc's VecGetArray acts differently in Fortran than it does in C.
911: ! Calling VecGetArray((Vec) X, (PetscScalar) x_array(0:1), (PetscOffset) x_index, info)
912: ! will return an array of doubles referenced by x_array offset by x_index.
913: ! i.e., to reference the kth element of X, use x_array(k + x_index).
914: ! Notice that by declaring the arrays with range (0:1), we are using the C 0-indexing practice.
915: PetscScalar xl_v(0:1)
916: PetscOffset xl_i
919: lb = -1.0d300
920: ub = 1.0d300
922: if (bmy .lt. 0) bmy = 0
923: if (bmy .gt. my) bmy = my
924: if (bmx .lt. 0) bmx = 0
925: if (bmx .gt. mx) bmx = mx
926:
928: call DAGetCorners(da,xs,ys,PETSC_NULL_INTEGER,xm,ym, &
929: & PETSC_NULL_INTEGER,info)
931: call VecSet(xl,lb,info)
932: call VecSet(xu,ub,info)
934: call VecGetArray(xl,xl_v,xl_i,info)
935:
937: do i=xs,xs+xm-1
939: do j=ys,ys+ym-1
940:
941: row=(j-ys)*xm + (i-xs)
943: if (i.ge.((mx-bmx)/2) .and. i.lt.(mx-(mx-bmx)/2) .and. &
944: & j.ge.((my-bmy)/2) .and. j.lt.(my-(my-bmy)/2)) then
945: xl_v(xl_i+row) = bheight
947: endif
949: enddo
950: enddo
953: call VecRestoreArray(xl,xl_v,xl_i,info)
954:
955: return
956: end
960:
961:
962: ! ----------------------------------------------------------------------------
963: !
964: !/*
965: ! MSA_InitialPoint - Calculates an obstacle for surface to stretch over
966: !
967: ! Output Parameter:
968: !. X - vector for initial guess
969: !
970: !*/
972: subroutine MSA_InitialPoint(X, info)
973: implicit none
975: ! mx,my,localX,da,Top,Left,Bottom,Right defined in plate2f.h
976: #include "plate2f.h"
977: Vec X
979: integer start,i,j,info
980: integer row,xs,xm,gxs,gxm,ys,ym,gys,gym
981: PetscScalar zero, np5
983: ! PETSc's VecGetArray acts differently in Fortran than it does in C.
984: ! Calling VecGetArray((Vec) X, (PetscScalar) x_array(0:1), (integer) x_index, info)
985: ! will return an array of doubles referenced by x_array offset by x_index.
986: ! i.e., to reference the kth element of X, use x_array(k + x_index).
987: ! Notice that by declaring the arrays with range (0:1), we are using the C 0-indexing practice.
988: PetscScalar left_v(0:1),right_v(0:1)
989: PetscScalar bottom_v(0:1),top_v(0:1)
990: PetscScalar x_v(0:1)
991: PetscOffset left_i, right_i, top_i
992: PetscOffset bottom_i,x_i
993: PetscTruth flg
994: PetscRandom rctx
995:
996: zero = 0.0d0
997: np5 = -0.5d0
999: call PetscOptionsGetInt(PETSC_NULL_CHARACTER,"-start", &
1000: & start,flg,info)
1002: if ((flg .eq. PETSC_TRUE) .and. (start .eq. 0)) then ! the zero vector is reasonable
1003: call VecSet(X,zero,info)
1005: elseif ((flg .eq. PETSC_TRUE) .and. (start .gt. 0)) then ! random start -0.5 < xi < 0.5
1006: call PetscRandomCreate(MPI_COMM_WORLD,rctx,info)
1007: do i=0,start-1
1008: call VecSetRandom(X,rctx,info)
1009: enddo
1011: call PetscRandomDestroy(rctx,info)
1012: call VecShift(X,np5,info)
1014: else ! average of boundary conditions
1015:
1016: ! Get Local mesh boundaries
1017: call DAGetCorners(da,xs,ys,PETSC_NULL_INTEGER,xm,ym, &
1018: & PETSC_NULL_INTEGER,info)
1019: call DAGetGhostCorners(da,gxs,gys,PETSC_NULL_INTEGER,gxm,gym, &
1020: & PETSC_NULL_INTEGER,info)
1024: ! Get pointers to vector data
1025: call VecGetArray(Top,top_v,top_i,info)
1026: call VecGetArray(Bottom,bottom_v,bottom_i,info)
1027: call VecGetArray(Left,left_v,left_i,info)
1028: call VecGetArray(Right,right_v,right_i,info)
1029:
1030: call VecGetArray(localX,x_v,x_i,info)
1031:
1032: ! Perform local computations
1033: do j=ys,ys+ym-1
1034: do i=xs,xs+xm-1
1035: row = (j-gys)*gxm + (i-gxs)
1036: x_v(x_i + row) = ((j+1)*bottom_v(bottom_i +i-xs+1)/my &
1037: & + (my-j+1)*top_v(top_i+i-xs+1)/(my+2) + &
1038: & (i+1)*left_v(left_i+j-ys+1)/mx + &
1039: & (mx-i+1)*right_v(right_i+j-ys+1)/(mx+2))*0.5
1040: enddo
1041: enddo
1043: ! Restore vectors
1044: call VecRestoreArray(localX,x_v,x_i,info)
1046: call VecRestoreArray(Left,left_v,left_i,info)
1047: call VecRestoreArray(Top,top_v,top_i,info)
1048: call VecRestoreArray(Bottom,bottom_v,bottom_i,info)
1049: call VecRestoreArray(Right,right_v,right_i,info)
1051: call DALocalToGlobal(da,localX,INSERT_VALUES,X,info)
1053: endif
1055: return
1056: end