Actual source code: ex48.c

petsc-dev 2014-04-23
Report Typos and Errors
  1: static const char help[] = "Toy hydrostatic ice flow with multigrid in 3D.\n\
  2: \n\
  3: Solves the hydrostatic (aka Blatter/Pattyn/First Order) equations for ice sheet flow\n\
  4: using multigrid.  The ice uses a power-law rheology with \"Glen\" exponent 3 (corresponds\n\
  5: to p=4/3 in a p-Laplacian).  The focus is on ISMIP-HOM experiments which assume periodic\n\
  6: boundary conditions in the x- and y-directions.\n\
  7: \n\
  8: Equations are rescaled so that the domain size and solution are O(1), details of this scaling\n\
  9: can be controlled by the options -units_meter, -units_second, and -units_kilogram.\n\
 10: \n\
 11: A VTK StructuredGrid output file can be written using the option -o filename.vts\n\
 12: \n\n";

 14: /*
 15: The equations for horizontal velocity (u,v) are

 17:   - [eta (4 u_x + 2 v_y)]_x - [eta (u_y + v_x)]_y - [eta u_z]_z + rho g s_x = 0
 18:   - [eta (4 v_y + 2 u_x)]_y - [eta (u_y + v_x)]_x - [eta v_z]_z + rho g s_y = 0

 20: where

 22:   eta = B/2 (epsilon + gamma)^((p-2)/2)

 24: is the nonlinear effective viscosity with regularization epsilon and hardness parameter B,
 25: written in terms of the second invariant

 27:   gamma = u_x^2 + v_y^2 + u_x v_y + (1/4) (u_y + v_x)^2 + (1/4) u_z^2 + (1/4) v_z^2

 29: The surface boundary conditions are the natural conditions.  The basal boundary conditions
 30: are either no-slip, or Navier (linear) slip with spatially variant friction coefficient beta^2.

 32: In the code, the equations for (u,v) are multiplied through by 1/(rho g) so that residuals are O(1).

 34: The discretization is Q1 finite elements, managed by a DMDA.  The grid is never distorted in the
 35: map (x,y) plane, but the bed and surface may be bumpy.  This is handled as usual in FEM, through
 36: the Jacobian of the coordinate transformation from a reference element to the physical element.

 38: Since ice-flow is tightly coupled in the z-direction (within columns), the DMDA is managed
 39: specially so that columns are never distributed, and are always contiguous in memory.
 40: This amounts to reversing the meaning of X,Y,Z compared to the DMDA's internal interpretation,
 41: and then indexing as vec[i][j][k].  The exotic coarse spaces require 2D DMDAs which are made to
 42: use compatible domain decomposition relative to the 3D DMDAs.

 44: There are two compile-time options:

 46:   NO_SSE2:
 47:     If the host supports SSE2, we use integration code that has been vectorized with SSE2
 48:     intrinsics, unless this macro is defined.  The intrinsics speed up integration by about
 49:     30% on my architecture (P8700, gcc-4.5 snapshot).

 51:   COMPUTE_LOWER_TRIANGULAR:
 52:     The element matrices we assemble are lower-triangular so it is not necessary to compute
 53:     all entries explicitly.  If this macro is defined, the lower-triangular entries are
 54:     computed explicitly.

 56: */

 58: #include <petscsnes.h>
 59: #include <ctype.h>              /* toupper() */
 60: #include <petsc-private/dmdaimpl.h>     /* There is not yet a public interface to manipulate dm->ops */

 62: #if !defined __STDC_VERSION__ || __STDC_VERSION__ < 199901L
 63: #  if defined __cplusplus       /* C++ restrict is nonstandard and compilers have inconsistent rules about where it can be used */
 64: #    define restrict
 65: #  else
 66: #    define restrict PETSC_RESTRICT
 67: #  endif
 68: #endif
 69: #if defined __SSE2__
 70: #  include <emmintrin.h>
 71: #endif

 73: /* The SSE2 kernels are only for PetscScalar=double on architectures that support it */
 74: #define USE_SSE2_KERNELS (!defined NO_SSE2                              \
 75:                           && !defined PETSC_USE_COMPLEX                 \
 76:                           && !defined PETSC_USE_REAL_SINGLE             \
 77:                           && !defined PETSC_USE_REAL___FLOAT128         \
 78:                           && defined __SSE2__)

 80: static PetscClassId THI_CLASSID;

 82: typedef enum {QUAD_GAUSS,QUAD_LOBATTO} QuadratureType;
 83: static const char      *QuadratureTypes[] = {"gauss","lobatto","QuadratureType","QUAD_",0};
 84: PETSC_UNUSED static const PetscReal HexQWeights[8]     = {1,1,1,1,1,1,1,1};
 85: PETSC_UNUSED static const PetscReal HexQNodes[]        = {-0.57735026918962573, 0.57735026918962573};
 86: #define G 0.57735026918962573
 87: #define H (0.5*(1.+G))
 88: #define L (0.5*(1.-G))
 89: #define M (-0.5)
 90: #define P (0.5)
 91: /* Special quadrature: Lobatto in horizontal, Gauss in vertical */
 92: static const PetscReal HexQInterp_Lobatto[8][8] = {{H,0,0,0,L,0,0,0},
 93:                                                    {0,H,0,0,0,L,0,0},
 94:                                                    {0,0,H,0,0,0,L,0},
 95:                                                    {0,0,0,H,0,0,0,L},
 96:                                                    {L,0,0,0,H,0,0,0},
 97:                                                    {0,L,0,0,0,H,0,0},
 98:                                                    {0,0,L,0,0,0,H,0},
 99:                                                    {0,0,0,L,0,0,0,H}};
100: static const PetscReal HexQDeriv_Lobatto[8][8][3] = {
101:   {{M*H,M*H,M},{P*H,0,0}  ,{0,0,0}    ,{0,P*H,0}  ,{M*L,M*L,P},{P*L,0,0}  ,{0,0,0}    ,{0,P*L,0}  },
102:   {{M*H,0,0}  ,{P*H,M*H,M},{0,P*H,0}  ,{0,0,0}    ,{M*L,0,0}  ,{P*L,M*L,P},{0,P*L,0}  ,{0,0,0}    },
103:   {{0,0,0}    ,{0,M*H,0}  ,{P*H,P*H,M},{M*H,0,0}  ,{0,0,0}    ,{0,M*L,0}  ,{P*L,P*L,P},{M*L,0,0}  },
104:   {{0,M*H,0}  ,{0,0,0}    ,{P*H,0,0}  ,{M*H,P*H,M},{0,M*L,0}  ,{0,0,0}    ,{P*L,0,0}  ,{M*L,P*L,P}},
105:   {{M*L,M*L,M},{P*L,0,0}  ,{0,0,0}    ,{0,P*L,0}  ,{M*H,M*H,P},{P*H,0,0}  ,{0,0,0}    ,{0,P*H,0}  },
106:   {{M*L,0,0}  ,{P*L,M*L,M},{0,P*L,0}  ,{0,0,0}    ,{M*H,0,0}  ,{P*H,M*H,P},{0,P*H,0}  ,{0,0,0}    },
107:   {{0,0,0}    ,{0,M*L,0}  ,{P*L,P*L,M},{M*L,0,0}  ,{0,0,0}    ,{0,M*H,0}  ,{P*H,P*H,P},{M*H,0,0}  },
108:   {{0,M*L,0}  ,{0,0,0}    ,{P*L,0,0}  ,{M*L,P*L,M},{0,M*H,0}  ,{0,0,0}    ,{P*H,0,0}  ,{M*H,P*H,P}}};
109: /* Stanndard Gauss */
110: static const PetscReal HexQInterp_Gauss[8][8] = {{H*H*H,L*H*H,L*L*H,H*L*H, H*H*L,L*H*L,L*L*L,H*L*L},
111:                                                  {L*H*H,H*H*H,H*L*H,L*L*H, L*H*L,H*H*L,H*L*L,L*L*L},
112:                                                  {L*L*H,H*L*H,H*H*H,L*H*H, L*L*L,H*L*L,H*H*L,L*H*L},
113:                                                  {H*L*H,L*L*H,L*H*H,H*H*H, H*L*L,L*L*L,L*H*L,H*H*L},
114:                                                  {H*H*L,L*H*L,L*L*L,H*L*L, H*H*H,L*H*H,L*L*H,H*L*H},
115:                                                  {L*H*L,H*H*L,H*L*L,L*L*L, L*H*H,H*H*H,H*L*H,L*L*H},
116:                                                  {L*L*L,H*L*L,H*H*L,L*H*L, L*L*H,H*L*H,H*H*H,L*H*H},
117:                                                  {H*L*L,L*L*L,L*H*L,H*H*L, H*L*H,L*L*H,L*H*H,H*H*H}};
118: static const PetscReal HexQDeriv_Gauss[8][8][3] = {
119:   {{M*H*H,H*M*H,H*H*M},{P*H*H,L*M*H,L*H*M},{P*L*H,L*P*H,L*L*M},{M*L*H,H*P*H,H*L*M}, {M*H*L,H*M*L,H*H*P},{P*H*L,L*M*L,L*H*P},{P*L*L,L*P*L,L*L*P},{M*L*L,H*P*L,H*L*P}},
120:   {{M*H*H,L*M*H,L*H*M},{P*H*H,H*M*H,H*H*M},{P*L*H,H*P*H,H*L*M},{M*L*H,L*P*H,L*L*M}, {M*H*L,L*M*L,L*H*P},{P*H*L,H*M*L,H*H*P},{P*L*L,H*P*L,H*L*P},{M*L*L,L*P*L,L*L*P}},
121:   {{M*L*H,L*M*H,L*L*M},{P*L*H,H*M*H,H*L*M},{P*H*H,H*P*H,H*H*M},{M*H*H,L*P*H,L*H*M}, {M*L*L,L*M*L,L*L*P},{P*L*L,H*M*L,H*L*P},{P*H*L,H*P*L,H*H*P},{M*H*L,L*P*L,L*H*P}},
122:   {{M*L*H,H*M*H,H*L*M},{P*L*H,L*M*H,L*L*M},{P*H*H,L*P*H,L*H*M},{M*H*H,H*P*H,H*H*M}, {M*L*L,H*M*L,H*L*P},{P*L*L,L*M*L,L*L*P},{P*H*L,L*P*L,L*H*P},{M*H*L,H*P*L,H*H*P}},
123:   {{M*H*L,H*M*L,H*H*M},{P*H*L,L*M*L,L*H*M},{P*L*L,L*P*L,L*L*M},{M*L*L,H*P*L,H*L*M}, {M*H*H,H*M*H,H*H*P},{P*H*H,L*M*H,L*H*P},{P*L*H,L*P*H,L*L*P},{M*L*H,H*P*H,H*L*P}},
124:   {{M*H*L,L*M*L,L*H*M},{P*H*L,H*M*L,H*H*M},{P*L*L,H*P*L,H*L*M},{M*L*L,L*P*L,L*L*M}, {M*H*H,L*M*H,L*H*P},{P*H*H,H*M*H,H*H*P},{P*L*H,H*P*H,H*L*P},{M*L*H,L*P*H,L*L*P}},
125:   {{M*L*L,L*M*L,L*L*M},{P*L*L,H*M*L,H*L*M},{P*H*L,H*P*L,H*H*M},{M*H*L,L*P*L,L*H*M}, {M*L*H,L*M*H,L*L*P},{P*L*H,H*M*H,H*L*P},{P*H*H,H*P*H,H*H*P},{M*H*H,L*P*H,L*H*P}},
126:   {{M*L*L,H*M*L,H*L*M},{P*L*L,L*M*L,L*L*M},{P*H*L,L*P*L,L*H*M},{M*H*L,H*P*L,H*H*M}, {M*L*H,H*M*H,H*L*P},{P*L*H,L*M*H,L*L*P},{P*H*H,L*P*H,L*H*P},{M*H*H,H*P*H,H*H*P}}};
127: static const PetscReal (*HexQInterp)[8],(*HexQDeriv)[8][3];
128: /* Standard 2x2 Gauss quadrature for the bottom layer. */
129: static const PetscReal QuadQInterp[4][4] = {{H*H,L*H,L*L,H*L},
130:                                             {L*H,H*H,H*L,L*L},
131:                                             {L*L,H*L,H*H,L*H},
132:                                             {H*L,L*L,L*H,H*H}};
133: static const PetscReal QuadQDeriv[4][4][2] = {
134:   {{M*H,M*H},{P*H,M*L},{P*L,P*L},{M*L,P*H}},
135:   {{M*H,M*L},{P*H,M*H},{P*L,P*H},{M*L,P*L}},
136:   {{M*L,M*L},{P*L,M*H},{P*H,P*H},{M*H,P*L}},
137:   {{M*L,M*H},{P*L,M*L},{P*H,P*L},{M*H,P*H}}};
138: #undef G
139: #undef H
140: #undef L
141: #undef M
142: #undef P

144: #define HexExtract(x,i,j,k,n) do {              \
145:     (n)[0] = (x)[i][j][k];                      \
146:     (n)[1] = (x)[i+1][j][k];                    \
147:     (n)[2] = (x)[i+1][j+1][k];                  \
148:     (n)[3] = (x)[i][j+1][k];                    \
149:     (n)[4] = (x)[i][j][k+1];                    \
150:     (n)[5] = (x)[i+1][j][k+1];                  \
151:     (n)[6] = (x)[i+1][j+1][k+1];                \
152:     (n)[7] = (x)[i][j+1][k+1];                  \
153:   } while (0)

155: #define HexExtractRef(x,i,j,k,n) do {           \
156:     (n)[0] = &(x)[i][j][k];                     \
157:     (n)[1] = &(x)[i+1][j][k];                   \
158:     (n)[2] = &(x)[i+1][j+1][k];                 \
159:     (n)[3] = &(x)[i][j+1][k];                   \
160:     (n)[4] = &(x)[i][j][k+1];                   \
161:     (n)[5] = &(x)[i+1][j][k+1];                 \
162:     (n)[6] = &(x)[i+1][j+1][k+1];               \
163:     (n)[7] = &(x)[i][j+1][k+1];                 \
164:   } while (0)

166: #define QuadExtract(x,i,j,n) do {               \
167:     (n)[0] = (x)[i][j];                         \
168:     (n)[1] = (x)[i+1][j];                       \
169:     (n)[2] = (x)[i+1][j+1];                     \
170:     (n)[3] = (x)[i][j+1];                       \
171:   } while (0)

173: static PetscScalar Sqr(PetscScalar a) {return a*a;}

175: static void HexGrad(const PetscReal dphi[][3],const PetscReal zn[],PetscReal dz[])
176: {
177:   PetscInt i;
178:   dz[0] = dz[1] = dz[2] = 0;
179:   for (i=0; i<8; i++) {
180:     dz[0] += dphi[i][0] * zn[i];
181:     dz[1] += dphi[i][1] * zn[i];
182:     dz[2] += dphi[i][2] * zn[i];
183:   }
184: }

186: static void HexComputeGeometry(PetscInt q,PetscReal hx,PetscReal hy,const PetscReal dz[restrict],PetscReal phi[restrict],PetscReal dphi[restrict][3],PetscReal *restrict jw)
187: {
188:   const PetscReal jac[3][3]  = {{hx/2,0,0}, {0,hy/2,0}, {dz[0],dz[1],dz[2]}};
189:   const PetscReal ijac[3][3] = {{1/jac[0][0],0,0}, {0,1/jac[1][1],0}, {-jac[2][0]/(jac[0][0]*jac[2][2]),-jac[2][1]/(jac[1][1]*jac[2][2]),1/jac[2][2]}};
190:   const PetscReal jdet       = jac[0][0]*jac[1][1]*jac[2][2];
191:   PetscInt        i;

193:   for (i=0; i<8; i++) {
194:     const PetscReal *dphir = HexQDeriv[q][i];
195:     phi[i]     = HexQInterp[q][i];
196:     dphi[i][0] = dphir[0]*ijac[0][0] + dphir[1]*ijac[1][0] + dphir[2]*ijac[2][0];
197:     dphi[i][1] = dphir[0]*ijac[0][1] + dphir[1]*ijac[1][1] + dphir[2]*ijac[2][1];
198:     dphi[i][2] = dphir[0]*ijac[0][2] + dphir[1]*ijac[1][2] + dphir[2]*ijac[2][2];
199:   }
200:   *jw = 1.0 * jdet;
201: }

203: typedef struct _p_THI   *THI;
204: typedef struct _n_Units *Units;

206: typedef struct {
207:   PetscScalar u,v;
208: } Node;

210: typedef struct {
211:   PetscScalar b;                /* bed */
212:   PetscScalar h;                /* thickness */
213:   PetscScalar beta2;            /* friction */
214: } PrmNode;

216: typedef struct {
217:   PetscReal min,max,cmin,cmax;
218: } PRange;

220: typedef enum {THIASSEMBLY_TRIDIAGONAL,THIASSEMBLY_FULL} THIAssemblyMode;

222: struct _p_THI {
223:   PETSCHEADER(int);
224:   void      (*initialize)(THI,PetscReal x,PetscReal y,PrmNode *p);
225:   PetscInt  zlevels;
226:   PetscReal Lx,Ly,Lz;           /* Model domain */
227:   PetscReal alpha;              /* Bed angle */
228:   Units     units;
229:   PetscReal dirichlet_scale;
230:   PetscReal ssa_friction_scale;
231:   PRange    eta;
232:   PRange    beta2;
233:   struct {
234:     PetscReal Bd2,eps,exponent;
235:   } viscosity;
236:   struct {
237:     PetscReal irefgam,eps2,exponent,refvel,epsvel;
238:   } friction;
239:   PetscReal rhog;
240:   PetscBool no_slip;
241:   PetscBool tridiagonal;
242:   PetscBool coarse2d;
243:   PetscBool verbose;
244:   MatType   mattype;
245: };

247: struct _n_Units {
248:   /* fundamental */
249:   PetscReal meter;
250:   PetscReal kilogram;
251:   PetscReal second;
252:   /* derived */
253:   PetscReal Pascal;
254:   PetscReal year;
255: };

257: static PetscErrorCode THIJacobianLocal_3D_Full(DMDALocalInfo*,Node***,Mat,Mat,THI);
258: static PetscErrorCode THIJacobianLocal_3D_Tridiagonal(DMDALocalInfo*,Node***,Mat,THI);
259: static PetscErrorCode THIJacobianLocal_2D(DMDALocalInfo*,Node**,Mat,THI);

261: static void PrmHexGetZ(const PrmNode pn[],PetscInt k,PetscInt zm,PetscReal zn[])
262: {
263:   const PetscScalar zm1    = zm-1,
264:                     znl[8] = {pn[0].b + pn[0].h*(PetscScalar)k/zm1,
265:                               pn[1].b + pn[1].h*(PetscScalar)k/zm1,
266:                               pn[2].b + pn[2].h*(PetscScalar)k/zm1,
267:                               pn[3].b + pn[3].h*(PetscScalar)k/zm1,
268:                               pn[0].b + pn[0].h*(PetscScalar)(k+1)/zm1,
269:                               pn[1].b + pn[1].h*(PetscScalar)(k+1)/zm1,
270:                               pn[2].b + pn[2].h*(PetscScalar)(k+1)/zm1,
271:                               pn[3].b + pn[3].h*(PetscScalar)(k+1)/zm1};
272:   PetscInt i;
273:   for (i=0; i<8; i++) zn[i] = PetscRealPart(znl[i]);
274: }

276: /* Tests A and C are from the ISMIP-HOM paper (Pattyn et al. 2008) */
277: static void THIInitialize_HOM_A(THI thi,PetscReal x,PetscReal y,PrmNode *p)
278: {
279:   Units     units = thi->units;
280:   PetscReal s     = -x*PetscSinReal(thi->alpha);

282:   p->b     = s - 1000*units->meter + 500*units->meter * PetscSinReal(x*2*PETSC_PI/thi->Lx) * PetscSinReal(y*2*PETSC_PI/thi->Ly);
283:   p->h     = s - p->b;
284:   p->beta2 = 1e30;
285: }

287: static void THIInitialize_HOM_C(THI thi,PetscReal x,PetscReal y,PrmNode *p)
288: {
289:   Units     units = thi->units;
290:   PetscReal s     = -x*PetscSinReal(thi->alpha);

292:   p->b = s - 1000*units->meter;
293:   p->h = s - p->b;
294:   /* tau_b = beta2 v   is a stress (Pa) */
295:   p->beta2 = 1000 * (1 + PetscSinReal(x*2*PETSC_PI/thi->Lx)*PetscSinReal(y*2*PETSC_PI/thi->Ly)) * units->Pascal * units->year / units->meter;
296: }

298: /* These are just toys */

300: /* Same bed as test A, free slip everywhere except for a discontinuous jump to a circular sticky region in the middle. */
301: static void THIInitialize_HOM_X(THI thi,PetscReal xx,PetscReal yy,PrmNode *p)
302: {
303:   Units     units = thi->units;
304:   PetscReal x     = xx*2*PETSC_PI/thi->Lx - PETSC_PI,y = yy*2*PETSC_PI/thi->Ly - PETSC_PI; /* [-pi,pi] */
305:   PetscReal r     = PetscSqrtReal(x*x + y*y),s = -x*PetscSinReal(thi->alpha);
306:   p->b     = s - 1000*units->meter + 500*units->meter*PetscSinReal(x + PETSC_PI) * PetscSinReal(y + PETSC_PI);
307:   p->h     = s - p->b;
308:   p->beta2 = 1000 * (r < 1 ? 2 : 0) * units->Pascal * units->year / units->meter;
309: }

311: /* Like Z, but with 200 meter cliffs */
312: static void THIInitialize_HOM_Y(THI thi,PetscReal xx,PetscReal yy,PrmNode *p)
313: {
314:   Units     units = thi->units;
315:   PetscReal x     = xx*2*PETSC_PI/thi->Lx - PETSC_PI,y = yy*2*PETSC_PI/thi->Ly - PETSC_PI; /* [-pi,pi] */
316:   PetscReal r     = PetscSqrtReal(x*x + y*y),s = -x*PetscSinReal(thi->alpha);

318:   p->b = s - 1000*units->meter + 500*units->meter * PetscSinReal(x + PETSC_PI) * PetscSinReal(y + PETSC_PI);
319:   if (PetscRealPart(p->b) > -700*units->meter) p->b += 200*units->meter;
320:   p->h     = s - p->b;
321:   p->beta2 = 1000 * (1. + PetscSinReal(PetscSqrtReal(16*r))/PetscSqrtReal(1e-2 + 16*r)*PetscCosReal(x*3/2)*PetscCosReal(y*3/2)) * units->Pascal * units->year / units->meter;
322: }

324: /* Same bed as A, smoothly varying slipperiness, similar to MATLAB's "sombrero" (uncorrelated with bathymetry) */
325: static void THIInitialize_HOM_Z(THI thi,PetscReal xx,PetscReal yy,PrmNode *p)
326: {
327:   Units     units = thi->units;
328:   PetscReal x     = xx*2*PETSC_PI/thi->Lx - PETSC_PI,y = yy*2*PETSC_PI/thi->Ly - PETSC_PI; /* [-pi,pi] */
329:   PetscReal r     = PetscSqrtReal(x*x + y*y),s = -x*PetscSinReal(thi->alpha);

331:   p->b     = s - 1000*units->meter + 500*units->meter * PetscSinReal(x + PETSC_PI) * PetscSinReal(y + PETSC_PI);
332:   p->h     = s - p->b;
333:   p->beta2 = 1000 * (1. + PetscSinReal(PetscSqrtReal(16*r))/PetscSqrtReal(1e-2 + 16*r)*PetscCosReal(x*3/2)*PetscCosReal(y*3/2)) * units->Pascal * units->year / units->meter;
334: }

336: static void THIFriction(THI thi,PetscReal rbeta2,PetscReal gam,PetscReal *beta2,PetscReal *dbeta2)
337: {
338:   if (thi->friction.irefgam == 0) {
339:     Units units = thi->units;
340:     thi->friction.irefgam = 1./(0.5*PetscSqr(thi->friction.refvel * units->meter / units->year));
341:     thi->friction.eps2    = 0.5*PetscSqr(thi->friction.epsvel * units->meter / units->year) * thi->friction.irefgam;
342:   }
343:   if (thi->friction.exponent == 0) {
344:     *beta2  = rbeta2;
345:     *dbeta2 = 0;
346:   } else {
347:     *beta2  = rbeta2 * PetscPowReal(thi->friction.eps2 + gam*thi->friction.irefgam,thi->friction.exponent);
348:     *dbeta2 = thi->friction.exponent * *beta2 / (thi->friction.eps2 + gam*thi->friction.irefgam) * thi->friction.irefgam;
349:   }
350: }

352: static void THIViscosity(THI thi,PetscReal gam,PetscReal *eta,PetscReal *deta)
353: {
354:   PetscReal Bd2,eps,exponent;
355:   if (thi->viscosity.Bd2 == 0) {
356:     Units units = thi->units;
357:     const PetscReal
358:       n = 3.,                                           /* Glen exponent */
359:       p = 1. + 1./n,                                    /* for Stokes */
360:       A = 1.e-16 * PetscPowReal(units->Pascal,-n) / units->year, /* softness parameter (Pa^{-n}/s) */
361:       B = PetscPowReal(A,-1./n);                                 /* hardness parameter */
362:     thi->viscosity.Bd2      = B/2;
363:     thi->viscosity.exponent = (p-2)/2;
364:     thi->viscosity.eps      = 0.5*PetscSqr(1e-5 / units->year);
365:   }
366:   Bd2      = thi->viscosity.Bd2;
367:   exponent = thi->viscosity.exponent;
368:   eps      = thi->viscosity.eps;
369:   *eta     = Bd2 * PetscPowReal(eps + gam,exponent);
370:   *deta    = exponent * (*eta) / (eps + gam);
371: }

373: static void RangeUpdate(PetscReal *min,PetscReal *max,PetscReal x)
374: {
375:   if (x < *min) *min = x;
376:   if (x > *max) *max = x;
377: }

379: static void PRangeClear(PRange *p)
380: {
381:   p->cmin = p->min = 1e100;
382:   p->cmax = p->max = -1e100;
383: }

387: static PetscErrorCode PRangeMinMax(PRange *p,PetscReal min,PetscReal max)
388: {

391:   p->cmin = min;
392:   p->cmax = max;
393:   if (min < p->min) p->min = min;
394:   if (max > p->max) p->max = max;
395:   return(0);
396: }

400: static PetscErrorCode THIDestroy(THI *thi)
401: {

405:   if (!*thi) return(0);
406:   if (--((PetscObject)(*thi))->refct > 0) {*thi = 0; return(0);}
407:   PetscFree((*thi)->units);
408:   PetscFree((*thi)->mattype);
409:   PetscHeaderDestroy(thi);
410:   return(0);
411: }

415: static PetscErrorCode THICreate(MPI_Comm comm,THI *inthi)
416: {
417:   static PetscBool registered = PETSC_FALSE;
418:   THI              thi;
419:   Units            units;
420:   PetscErrorCode   ierr;

423:   *inthi = 0;
424:   if (!registered) {
425:     PetscClassIdRegister("Toy Hydrostatic Ice",&THI_CLASSID);
426:     registered = PETSC_TRUE;
427:   }
428:   PetscHeaderCreate(thi,_p_THI,0,THI_CLASSID,"THI","Toy Hydrostatic Ice","",comm,THIDestroy,0);

430:   PetscNew(&thi->units);
431:   units           = thi->units;
432:   units->meter    = 1e-2;
433:   units->second   = 1e-7;
434:   units->kilogram = 1e-12;

436:   PetscOptionsBegin(comm,NULL,"Scaled units options","");
437:   {
438:     PetscOptionsReal("-units_meter","1 meter in scaled length units","",units->meter,&units->meter,NULL);
439:     PetscOptionsReal("-units_second","1 second in scaled time units","",units->second,&units->second,NULL);
440:     PetscOptionsReal("-units_kilogram","1 kilogram in scaled mass units","",units->kilogram,&units->kilogram,NULL);
441:   }
442:   PetscOptionsEnd();
443:   units->Pascal = units->kilogram / (units->meter * PetscSqr(units->second));
444:   units->year   = 31556926. * units->second, /* seconds per year */

446:   thi->Lx              = 10.e3;
447:   thi->Ly              = 10.e3;
448:   thi->Lz              = 1000;
449:   thi->dirichlet_scale = 1;
450:   thi->verbose         = PETSC_FALSE;

452:   PetscOptionsBegin(comm,NULL,"Toy Hydrostatic Ice options","");
453:   {
454:     QuadratureType quad       = QUAD_GAUSS;
455:     char           homexp[]   = "A";
456:     char           mtype[256] = MATSBAIJ;
457:     PetscReal      L,m = 1.0;
458:     PetscBool      flg;
459:     L    = thi->Lx;
460:     PetscOptionsReal("-thi_L","Domain size (m)","",L,&L,&flg);
461:     if (flg) thi->Lx = thi->Ly = L;
462:     PetscOptionsReal("-thi_Lx","X Domain size (m)","",thi->Lx,&thi->Lx,NULL);
463:     PetscOptionsReal("-thi_Ly","Y Domain size (m)","",thi->Ly,&thi->Ly,NULL);
464:     PetscOptionsReal("-thi_Lz","Z Domain size (m)","",thi->Lz,&thi->Lz,NULL);
465:     PetscOptionsString("-thi_hom","ISMIP-HOM experiment (A or C)","",homexp,homexp,sizeof(homexp),NULL);
466:     switch (homexp[0] = toupper(homexp[0])) {
467:     case 'A':
468:       thi->initialize = THIInitialize_HOM_A;
469:       thi->no_slip    = PETSC_TRUE;
470:       thi->alpha      = 0.5;
471:       break;
472:     case 'C':
473:       thi->initialize = THIInitialize_HOM_C;
474:       thi->no_slip    = PETSC_FALSE;
475:       thi->alpha      = 0.1;
476:       break;
477:     case 'X':
478:       thi->initialize = THIInitialize_HOM_X;
479:       thi->no_slip    = PETSC_FALSE;
480:       thi->alpha      = 0.3;
481:       break;
482:     case 'Y':
483:       thi->initialize = THIInitialize_HOM_Y;
484:       thi->no_slip    = PETSC_FALSE;
485:       thi->alpha      = 0.5;
486:       break;
487:     case 'Z':
488:       thi->initialize = THIInitialize_HOM_Z;
489:       thi->no_slip    = PETSC_FALSE;
490:       thi->alpha      = 0.5;
491:       break;
492:     default:
493:       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"HOM experiment '%c' not implemented",homexp[0]);
494:     }
495:     PetscOptionsEnum("-thi_quadrature","Quadrature to use for 3D elements","",QuadratureTypes,(PetscEnum)quad,(PetscEnum*)&quad,NULL);
496:     switch (quad) {
497:     case QUAD_GAUSS:
498:       HexQInterp = HexQInterp_Gauss;
499:       HexQDeriv  = HexQDeriv_Gauss;
500:       break;
501:     case QUAD_LOBATTO:
502:       HexQInterp = HexQInterp_Lobatto;
503:       HexQDeriv  = HexQDeriv_Lobatto;
504:       break;
505:     }
506:     PetscOptionsReal("-thi_alpha","Bed angle (degrees)","",thi->alpha,&thi->alpha,NULL);

508:     thi->friction.refvel = 100.;
509:     thi->friction.epsvel = 1.;

511:     PetscOptionsReal("-thi_friction_refvel","Reference velocity for sliding","",thi->friction.refvel,&thi->friction.refvel,NULL);
512:     PetscOptionsReal("-thi_friction_epsvel","Regularization velocity for sliding","",thi->friction.epsvel,&thi->friction.epsvel,NULL);
513:     PetscOptionsReal("-thi_friction_m","Friction exponent, 0=Coulomb, 1=Navier","",m,&m,NULL);

515:     thi->friction.exponent = (m-1)/2;

517:     PetscOptionsReal("-thi_dirichlet_scale","Scale Dirichlet boundary conditions by this factor","",thi->dirichlet_scale,&thi->dirichlet_scale,NULL);
518:     PetscOptionsReal("-thi_ssa_friction_scale","Scale slip boundary conditions by this factor in SSA (2D) assembly","",thi->ssa_friction_scale,&thi->ssa_friction_scale,NULL);
519:     PetscOptionsBool("-thi_coarse2d","Use a 2D coarse space corresponding to SSA","",thi->coarse2d,&thi->coarse2d,NULL);
520:     PetscOptionsBool("-thi_tridiagonal","Assemble a tridiagonal system (column coupling only) on the finest level","",thi->tridiagonal,&thi->tridiagonal,NULL);
521:     PetscOptionsFList("-thi_mat_type","Matrix type","MatSetType",MatList,mtype,(char*)mtype,sizeof(mtype),NULL);
522:     PetscStrallocpy(mtype,(char**)&thi->mattype);
523:     PetscOptionsBool("-thi_verbose","Enable verbose output (like matrix sizes and statistics)","",thi->verbose,&thi->verbose,NULL);
524:   }
525:   PetscOptionsEnd();

527:   /* dimensionalize */
528:   thi->Lx    *= units->meter;
529:   thi->Ly    *= units->meter;
530:   thi->Lz    *= units->meter;
531:   thi->alpha *= PETSC_PI / 180;

533:   PRangeClear(&thi->eta);
534:   PRangeClear(&thi->beta2);

536:   {
537:     PetscReal u       = 1000*units->meter/(3e7*units->second),
538:               gradu   = u / (100*units->meter),eta,deta,
539:               rho     = 910 * units->kilogram/PetscPowReal(units->meter,3),
540:               grav    = 9.81 * units->meter/PetscSqr(units->second),
541:               driving = rho * grav * PetscSinReal(thi->alpha) * 1000*units->meter;
542:     THIViscosity(thi,0.5*gradu*gradu,&eta,&deta);
543:     thi->rhog = rho * grav;
544:     if (thi->verbose) {
545:       PetscPrintf(PetscObjectComm((PetscObject)thi),"Units: meter %8.2g  second %8.2g  kg %8.2g  Pa %8.2g\n",(double)units->meter,(double)units->second,(double)units->kilogram,(double)units->Pascal);
546:       PetscPrintf(PetscObjectComm((PetscObject)thi),"Domain (%6.2g,%6.2g,%6.2g), pressure %8.2g, driving stress %8.2g\n",(double)thi->Lx,(double)thi->Ly,(double)thi->Lz,(double)(rho*grav*1e3*units->meter),(double)driving);
547:       PetscPrintf(PetscObjectComm((PetscObject)thi),"Large velocity 1km/a %8.2g, velocity gradient %8.2g, eta %8.2g, stress %8.2g, ratio %8.2g\n",(double)u,(double)gradu,(double)eta,(double)(2*eta*gradu),(double)(2*eta*gradu/driving));
548:       THIViscosity(thi,0.5*PetscSqr(1e-3*gradu),&eta,&deta);
549:       PetscPrintf(PetscObjectComm((PetscObject)thi),"Small velocity 1m/a  %8.2g, velocity gradient %8.2g, eta %8.2g, stress %8.2g, ratio %8.2g\n",(double)(1e-3*u),(double)(1e-3*gradu),(double)eta,(double)(2*eta*1e-3*gradu),(double)(2*eta*1e-3*gradu/driving));
550:     }
551:   }

553:   *inthi = thi;
554:   return(0);
555: }

559: static PetscErrorCode THIInitializePrm(THI thi,DM da2prm,Vec prm)
560: {
561:   PrmNode        **p;
562:   PetscInt       i,j,xs,xm,ys,ym,mx,my;

566:   DMDAGetGhostCorners(da2prm,&ys,&xs,0,&ym,&xm,0);
567:   DMDAGetInfo(da2prm,0, &my,&mx,0, 0,0,0, 0,0,0,0,0,0);
568:   DMDAVecGetArray(da2prm,prm,&p);
569:   for (i=xs; i<xs+xm; i++) {
570:     for (j=ys; j<ys+ym; j++) {
571:       PetscReal xx = thi->Lx*i/mx,yy = thi->Ly*j/my;
572:       thi->initialize(thi,xx,yy,&p[i][j]);
573:     }
574:   }
575:   DMDAVecRestoreArray(da2prm,prm,&p);
576:   return(0);
577: }

581: static PetscErrorCode THISetUpDM(THI thi,DM dm)
582: {
583:   PetscErrorCode  ierr;
584:   PetscInt        refinelevel,coarsenlevel,level,dim,Mx,My,Mz,mx,my,s;
585:   DMDAStencilType st;
586:   DM              da2prm;
587:   Vec             X;

590:   DMDAGetInfo(dm,&dim, &Mz,&My,&Mx, 0,&my,&mx, 0,&s,0,0,0,&st);
591:   if (dim == 2) {
592:     DMDAGetInfo(dm,&dim, &My,&Mx,0, &my,&mx,0, 0,&s,0,0,0,&st);
593:   }
594:   DMGetRefineLevel(dm,&refinelevel);
595:   DMGetCoarsenLevel(dm,&coarsenlevel);
596:   level = refinelevel - coarsenlevel;
597:   DMDACreate2d(PetscObjectComm((PetscObject)thi),DM_BOUNDARY_PERIODIC,DM_BOUNDARY_PERIODIC,st,My,Mx,my,mx,sizeof(PrmNode)/sizeof(PetscScalar),s,0,0,&da2prm);
598:   DMCreateLocalVector(da2prm,&X);
599:   {
600:     PetscReal Lx = thi->Lx / thi->units->meter,Ly = thi->Ly / thi->units->meter,Lz = thi->Lz / thi->units->meter;
601:     if (dim == 2) {
602:       PetscPrintf(PetscObjectComm((PetscObject)thi),"Level %D domain size (m) %8.2g x %8.2g, num elements %D x %D (%D), size (m) %g x %g\n",level,(double)Lx,(double)Ly,Mx,My,Mx*My,(double)(Lx/Mx),(double)(Ly/My));
603:     } else {
604:       PetscPrintf(PetscObjectComm((PetscObject)thi),"Level %D domain size (m) %8.2g x %8.2g x %8.2g, num elements %D x %D x %D (%D), size (m) %g x %g x %g\n",level,(double)Lx,(double)Ly,(double)Lz,Mx,My,Mz,Mx*My*Mz,(double)(Lx/Mx),(double)(Ly/My),(double)(1000./(Mz-1)));
605:     }
606:   }
607:   THIInitializePrm(thi,da2prm,X);
608:   if (thi->tridiagonal) {       /* Reset coarse Jacobian evaluation */
609:     DMDASNESSetJacobianLocal(dm,(DMDASNESJacobian)THIJacobianLocal_3D_Full,thi);
610:   }
611:   if (thi->coarse2d) {
612:     DMDASNESSetJacobianLocal(dm,(DMDASNESJacobian)THIJacobianLocal_2D,thi);
613:   }
614:   PetscObjectCompose((PetscObject)dm,"DMDA2Prm",(PetscObject)da2prm);
615:   PetscObjectCompose((PetscObject)dm,"DMDA2Prm_Vec",(PetscObject)X);
616:   DMDestroy(&da2prm);
617:   VecDestroy(&X);
618:   return(0);
619: }

623: static PetscErrorCode DMCoarsenHook_THI(DM dmf,DM dmc,void *ctx)
624: {
625:   THI            thi = (THI)ctx;
627:   PetscInt       rlevel,clevel;

630:   THISetUpDM(thi,dmc);
631:   DMGetRefineLevel(dmc,&rlevel);
632:   DMGetCoarsenLevel(dmc,&clevel);
633:   if (rlevel-clevel == 0) {DMSetMatType(dmc,MATAIJ);}
634:   DMCoarsenHookAdd(dmc,DMCoarsenHook_THI,NULL,thi);
635:   return(0);
636: }

640: static PetscErrorCode DMRefineHook_THI(DM dmc,DM dmf,void *ctx)
641: {
642:   THI            thi = (THI)ctx;

646:   THISetUpDM(thi,dmf);
647:   DMSetMatType(dmf,thi->mattype);
648:   DMRefineHookAdd(dmf,DMRefineHook_THI,NULL,thi);
649:   /* With grid sequencing, a formerly-refined DM will later be coarsened by PCSetUp_MG */
650:   DMCoarsenHookAdd(dmf,DMCoarsenHook_THI,NULL,thi);
651:   return(0);
652: }

656: static PetscErrorCode THIDAGetPrm(DM da,PrmNode ***prm)
657: {
659:   DM             da2prm;
660:   Vec            X;

663:   PetscObjectQuery((PetscObject)da,"DMDA2Prm",(PetscObject*)&da2prm);
664:   if (!da2prm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"No DMDA2Prm composed with given DMDA");
665:   PetscObjectQuery((PetscObject)da,"DMDA2Prm_Vec",(PetscObject*)&X);
666:   if (!X) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"No DMDA2Prm_Vec composed with given DMDA");
667:   DMDAVecGetArray(da2prm,X,prm);
668:   return(0);
669: }

673: static PetscErrorCode THIDARestorePrm(DM da,PrmNode ***prm)
674: {
676:   DM             da2prm;
677:   Vec            X;

680:   PetscObjectQuery((PetscObject)da,"DMDA2Prm",(PetscObject*)&da2prm);
681:   if (!da2prm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"No DMDA2Prm composed with given DMDA");
682:   PetscObjectQuery((PetscObject)da,"DMDA2Prm_Vec",(PetscObject*)&X);
683:   if (!X) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"No DMDA2Prm_Vec composed with given DMDA");
684:   DMDAVecRestoreArray(da2prm,X,prm);
685:   return(0);
686: }

690: static PetscErrorCode THIInitial(SNES snes,Vec X,void *ctx)
691: {
692:   THI            thi;
693:   PetscInt       i,j,k,xs,xm,ys,ym,zs,zm,mx,my;
694:   PetscReal      hx,hy;
695:   PrmNode        **prm;
696:   Node           ***x;
698:   DM             da;

701:   SNESGetDM(snes,&da);
702:   DMGetApplicationContext(da,&thi);
703:   DMDAGetInfo(da,0, 0,&my,&mx, 0,0,0, 0,0,0,0,0,0);
704:   DMDAGetCorners(da,&zs,&ys,&xs,&zm,&ym,&xm);
705:   DMDAVecGetArray(da,X,&x);
706:   THIDAGetPrm(da,&prm);
707:   hx   = thi->Lx / mx;
708:   hy   = thi->Ly / my;
709:   for (i=xs; i<xs+xm; i++) {
710:     for (j=ys; j<ys+ym; j++) {
711:       for (k=zs; k<zs+zm; k++) {
712:         const PetscScalar zm1      = zm-1,
713:                           drivingx = thi->rhog * (prm[i+1][j].b+prm[i+1][j].h - prm[i-1][j].b-prm[i-1][j].h) / (2*hx),
714:                           drivingy = thi->rhog * (prm[i][j+1].b+prm[i][j+1].h - prm[i][j-1].b-prm[i][j-1].h) / (2*hy);
715:         x[i][j][k].u = 0. * drivingx * prm[i][j].h*(PetscScalar)k/zm1;
716:         x[i][j][k].v = 0. * drivingy * prm[i][j].h*(PetscScalar)k/zm1;
717:       }
718:     }
719:   }
720:   DMDAVecRestoreArray(da,X,&x);
721:   THIDARestorePrm(da,&prm);
722:   return(0);
723: }

725: static void PointwiseNonlinearity(THI thi,const Node n[restrict],const PetscReal phi[restrict],PetscReal dphi[restrict][3],PetscScalar *restrict u,PetscScalar *restrict v,PetscScalar du[restrict],PetscScalar dv[restrict],PetscReal *eta,PetscReal *deta)
726: {
727:   PetscInt    l,ll;
728:   PetscScalar gam;

730:   du[0] = du[1] = du[2] = 0;
731:   dv[0] = dv[1] = dv[2] = 0;
732:   *u    = 0;
733:   *v    = 0;
734:   for (l=0; l<8; l++) {
735:     *u += phi[l] * n[l].u;
736:     *v += phi[l] * n[l].v;
737:     for (ll=0; ll<3; ll++) {
738:       du[ll] += dphi[l][ll] * n[l].u;
739:       dv[ll] += dphi[l][ll] * n[l].v;
740:     }
741:   }
742:   gam = Sqr(du[0]) + Sqr(dv[1]) + du[0]*dv[1] + 0.25*Sqr(du[1]+dv[0]) + 0.25*Sqr(du[2]) + 0.25*Sqr(dv[2]);
743:   THIViscosity(thi,PetscRealPart(gam),eta,deta);
744: }

746: static void PointwiseNonlinearity2D(THI thi,Node n[],PetscReal phi[],PetscReal dphi[4][2],PetscScalar *u,PetscScalar *v,PetscScalar du[],PetscScalar dv[],PetscReal *eta,PetscReal *deta)
747: {
748:   PetscInt    l,ll;
749:   PetscScalar gam;

751:   du[0] = du[1] = 0;
752:   dv[0] = dv[1] = 0;
753:   *u    = 0;
754:   *v    = 0;
755:   for (l=0; l<4; l++) {
756:     *u += phi[l] * n[l].u;
757:     *v += phi[l] * n[l].v;
758:     for (ll=0; ll<2; ll++) {
759:       du[ll] += dphi[l][ll] * n[l].u;
760:       dv[ll] += dphi[l][ll] * n[l].v;
761:     }
762:   }
763:   gam = Sqr(du[0]) + Sqr(dv[1]) + du[0]*dv[1] + 0.25*Sqr(du[1]+dv[0]);
764:   THIViscosity(thi,PetscRealPart(gam),eta,deta);
765: }

769: static PetscErrorCode THIFunctionLocal(DMDALocalInfo *info,Node ***x,Node ***f,THI thi)
770: {
771:   PetscInt       xs,ys,xm,ym,zm,i,j,k,q,l;
772:   PetscReal      hx,hy,etamin,etamax,beta2min,beta2max;
773:   PrmNode        **prm;

777:   xs = info->zs;
778:   ys = info->ys;
779:   xm = info->zm;
780:   ym = info->ym;
781:   zm = info->xm;
782:   hx = thi->Lx / info->mz;
783:   hy = thi->Ly / info->my;

785:   etamin   = 1e100;
786:   etamax   = 0;
787:   beta2min = 1e100;
788:   beta2max = 0;

790:   THIDAGetPrm(info->da,&prm);

792:   for (i=xs; i<xs+xm; i++) {
793:     for (j=ys; j<ys+ym; j++) {
794:       PrmNode pn[4];
795:       QuadExtract(prm,i,j,pn);
796:       for (k=0; k<zm-1; k++) {
797:         PetscInt  ls = 0;
798:         Node      n[8],*fn[8];
799:         PetscReal zn[8],etabase = 0;
800:         PrmHexGetZ(pn,k,zm,zn);
801:         HexExtract(x,i,j,k,n);
802:         HexExtractRef(f,i,j,k,fn);
803:         if (thi->no_slip && k == 0) {
804:           for (l=0; l<4; l++) n[l].u = n[l].v = 0;
805:           /* The first 4 basis functions lie on the bottom layer, so their contribution is exactly 0, hence we can skip them */
806:           ls = 4;
807:         }
808:         for (q=0; q<8; q++) {
809:           PetscReal   dz[3],phi[8],dphi[8][3],jw,eta,deta;
810:           PetscScalar du[3],dv[3],u,v;
811:           HexGrad(HexQDeriv[q],zn,dz);
812:           HexComputeGeometry(q,hx,hy,dz,phi,dphi,&jw);
813:           PointwiseNonlinearity(thi,n,phi,dphi,&u,&v,du,dv,&eta,&deta);
814:           jw /= thi->rhog;      /* scales residuals to be O(1) */
815:           if (q == 0) etabase = eta;
816:           RangeUpdate(&etamin,&etamax,eta);
817:           for (l=ls; l<8; l++) { /* test functions */
818:             const PetscReal ds[2] = {-PetscSinReal(thi->alpha),0};
819:             const PetscReal pp    = phi[l],*dp = dphi[l];
820:             fn[l]->u += dp[0]*jw*eta*(4.*du[0]+2.*dv[1]) + dp[1]*jw*eta*(du[1]+dv[0]) + dp[2]*jw*eta*du[2] + pp*jw*thi->rhog*ds[0];
821:             fn[l]->v += dp[1]*jw*eta*(2.*du[0]+4.*dv[1]) + dp[0]*jw*eta*(du[1]+dv[0]) + dp[2]*jw*eta*dv[2] + pp*jw*thi->rhog*ds[1];
822:           }
823:         }
824:         if (k == 0) { /* we are on a bottom face */
825:           if (thi->no_slip) {
826:             /* Note: Non-Galerkin coarse grid operators are very sensitive to the scaling of Dirichlet boundary
827:             * conditions.  After shenanigans above, etabase contains the effective viscosity at the closest quadrature
828:             * point to the bed.  We want the diagonal entry in the Dirichlet condition to have similar magnitude to the
829:             * diagonal entry corresponding to the adjacent node.  The fundamental scaling of the viscous part is in
830:             * diagu, diagv below.  This scaling is easy to recognize by considering the finite difference operator after
831:             * scaling by element size.  The no-slip Dirichlet condition is scaled by this factor, and also in the
832:             * assembled matrix (see the similar block in THIJacobianLocal).
833:             *
834:             * Note that the residual at this Dirichlet node is linear in the state at this node, but also depends
835:             * (nonlinearly in general) on the neighboring interior nodes through the local viscosity.  This will make
836:             * a matrix-free Jacobian have extra entries in the corresponding row.  We assemble only the diagonal part,
837:             * so the solution will exactly satisfy the boundary condition after the first linear iteration.
838:             */
839:             const PetscReal   hz    = PetscRealPart(pn[0].h)/(zm-1.);
840:             const PetscScalar diagu = 2*etabase/thi->rhog*(hx*hy/hz + hx*hz/hy + 4*hy*hz/hx),diagv = 2*etabase/thi->rhog*(hx*hy/hz + 4*hx*hz/hy + hy*hz/hx);
841:             fn[0]->u = thi->dirichlet_scale*diagu*x[i][j][k].u;
842:             fn[0]->v = thi->dirichlet_scale*diagv*x[i][j][k].v;
843:           } else {              /* Integrate over bottom face to apply boundary condition */
844:             for (q=0; q<4; q++) {
845:               const PetscReal jw = 0.25*hx*hy/thi->rhog,*phi = QuadQInterp[q];
846:               PetscScalar     u  =0,v=0,rbeta2=0;
847:               PetscReal       beta2,dbeta2;
848:               for (l=0; l<4; l++) {
849:                 u      += phi[l]*n[l].u;
850:                 v      += phi[l]*n[l].v;
851:                 rbeta2 += phi[l]*pn[l].beta2;
852:               }
853:               THIFriction(thi,PetscRealPart(rbeta2),PetscRealPart(u*u+v*v)/2,&beta2,&dbeta2);
854:               RangeUpdate(&beta2min,&beta2max,beta2);
855:               for (l=0; l<4; l++) {
856:                 const PetscReal pp = phi[l];
857:                 fn[ls+l]->u += pp*jw*beta2*u;
858:                 fn[ls+l]->v += pp*jw*beta2*v;
859:               }
860:             }
861:           }
862:         }
863:       }
864:     }
865:   }

867:   THIDARestorePrm(info->da,&prm);

869:   PRangeMinMax(&thi->eta,etamin,etamax);
870:   PRangeMinMax(&thi->beta2,beta2min,beta2max);
871:   return(0);
872: }

876: static PetscErrorCode THIMatrixStatistics(THI thi,Mat B,PetscViewer viewer)
877: {
879:   PetscReal      nrm;
880:   PetscInt       m;
881:   PetscMPIInt    rank;

884:   MatNorm(B,NORM_FROBENIUS,&nrm);
885:   MatGetSize(B,&m,0);
886:   MPI_Comm_rank(PetscObjectComm((PetscObject)B),&rank);
887:   if (!rank) {
888:     PetscScalar val0,val2;
889:     MatGetValue(B,0,0,&val0);
890:     MatGetValue(B,2,2,&val2);
891:     PetscViewerASCIIPrintf(viewer,"Matrix dim %D norm %8.2e (0,0) %8.2e  (2,2) %8.2e %8.2e <= eta <= %8.2e %8.2e <= beta2 <= %8.2e\n",m,(double)nrm,(double)PetscRealPart(val0),(double)PetscRealPart(val2),(double)thi->eta.cmin,(double)thi->eta.cmax,(double)thi->beta2.cmin,(double)thi->beta2.cmax);
892:   }
893:   return(0);
894: }

898: static PetscErrorCode THISurfaceStatistics(DM da,Vec X,PetscReal *min,PetscReal *max,PetscReal *mean)
899: {
901:   Node           ***x;
902:   PetscInt       i,j,xs,ys,zs,xm,ym,zm,mx,my,mz;
903:   PetscReal      umin = 1e100,umax=-1e100;
904:   PetscScalar    usum = 0.0,gusum;

907:   *min = *max = *mean = 0;
908:   DMDAGetInfo(da,0, &mz,&my,&mx, 0,0,0, 0,0,0,0,0,0);
909:   DMDAGetCorners(da,&zs,&ys,&xs,&zm,&ym,&xm);
910:   if (zs != 0 || zm != mz) SETERRQ(PETSC_COMM_SELF,1,"Unexpected decomposition");
911:   DMDAVecGetArray(da,X,&x);
912:   for (i=xs; i<xs+xm; i++) {
913:     for (j=ys; j<ys+ym; j++) {
914:       PetscReal u = PetscRealPart(x[i][j][zm-1].u);
915:       RangeUpdate(&umin,&umax,u);
916:       usum += u;
917:     }
918:   }
919:   DMDAVecRestoreArray(da,X,&x);
920:   MPI_Allreduce(&umin,min,1,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)da));
921:   MPI_Allreduce(&umax,max,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)da));
922:   MPI_Allreduce(&usum,&gusum,1,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)da));
923:   *mean = PetscRealPart(gusum) / (mx*my);
924:   return(0);
925: }

929: static PetscErrorCode THISolveStatistics(THI thi,SNES snes,PetscInt coarsened,const char name[])
930: {
931:   MPI_Comm       comm;
932:   Vec            X;
933:   DM             dm;

937:   PetscObjectGetComm((PetscObject)thi,&comm);
938:   SNESGetSolution(snes,&X);
939:   SNESGetDM(snes,&dm);
940:   PetscPrintf(comm,"Solution statistics after solve: %s\n",name);
941:   {
942:     PetscInt            its,lits;
943:     SNESConvergedReason reason;
944:     SNESGetIterationNumber(snes,&its);
945:     SNESGetConvergedReason(snes,&reason);
946:     SNESGetLinearSolveIterations(snes,&lits);
947:     PetscPrintf(comm,"%s: Number of SNES iterations = %D, total linear iterations = %D\n",SNESConvergedReasons[reason],its,lits);
948:   }
949:   {
950:     PetscReal   nrm2,tmin[3]={1e100,1e100,1e100},tmax[3]={-1e100,-1e100,-1e100},min[3],max[3];
951:     PetscInt    i,j,m;
952:     PetscScalar *x;
953:     VecNorm(X,NORM_2,&nrm2);
954:     VecGetLocalSize(X,&m);
955:     VecGetArray(X,&x);
956:     for (i=0; i<m; i+=2) {
957:       PetscReal u = PetscRealPart(x[i]),v = PetscRealPart(x[i+1]),c = PetscSqrtReal(u*u+v*v);
958:       tmin[0] = PetscMin(u,tmin[0]);
959:       tmin[1] = PetscMin(v,tmin[1]);
960:       tmin[2] = PetscMin(c,tmin[2]);
961:       tmax[0] = PetscMax(u,tmax[0]);
962:       tmax[1] = PetscMax(v,tmax[1]);
963:       tmax[2] = PetscMax(c,tmax[2]);
964:     }
965:     VecRestoreArray(X,&x);
966:     MPI_Allreduce(tmin,min,3,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)thi));
967:     MPI_Allreduce(tmax,max,3,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)thi));
968:     /* Dimensionalize to meters/year */
969:     nrm2 *= thi->units->year / thi->units->meter;
970:     for (j=0; j<3; j++) {
971:       min[j] *= thi->units->year / thi->units->meter;
972:       max[j] *= thi->units->year / thi->units->meter;
973:     }
974:     PetscPrintf(comm,"|X|_2 %g   %g <= u <=  %g   %g <= v <=  %g   %g <= c <=  %g \n",(double)nrm2,(double)min[0],(double)max[0],(double)min[1],(double)max[1],(double)min[2],(double)max[2]);
975:     {
976:       PetscReal umin,umax,umean;
977:       THISurfaceStatistics(dm,X,&umin,&umax,&umean);
978:       umin  *= thi->units->year / thi->units->meter;
979:       umax  *= thi->units->year / thi->units->meter;
980:       umean *= thi->units->year / thi->units->meter;
981:       PetscPrintf(comm,"Surface statistics: u in [%12.6e, %12.6e] mean %12.6e\n",(double)umin,(double)umax,(double)umean);
982:     }
983:     /* These values stay nondimensional */
984:     PetscPrintf(comm,"Global eta range   %g to %g converged range %g to %g\n",(double)thi->eta.min,(double)thi->eta.max,(double)thi->eta.cmin,(double)thi->eta.cmax);
985:     PetscPrintf(comm,"Global beta2 range %g to %g converged range %g to %g\n",(double)thi->beta2.min,(double)thi->beta2.max,(double)thi->beta2.cmin,(double)thi->beta2.cmax);
986:   }
987:   return(0);
988: }

992: static PetscErrorCode THIJacobianLocal_2D(DMDALocalInfo *info,Node **x,Mat B,THI thi)
993: {
994:   PetscInt       xs,ys,xm,ym,i,j,q,l,ll;
995:   PetscReal      hx,hy;
996:   PrmNode        **prm;

1000:   xs = info->ys;
1001:   ys = info->xs;
1002:   xm = info->ym;
1003:   ym = info->xm;
1004:   hx = thi->Lx / info->my;
1005:   hy = thi->Ly / info->mx;

1007:   MatZeroEntries(B);
1008:   THIDAGetPrm(info->da,&prm);

1010:   for (i=xs; i<xs+xm; i++) {
1011:     for (j=ys; j<ys+ym; j++) {
1012:       Node        n[4];
1013:       PrmNode     pn[4];
1014:       PetscScalar Ke[4*2][4*2];
1015:       QuadExtract(prm,i,j,pn);
1016:       QuadExtract(x,i,j,n);
1017:       PetscMemzero(Ke,sizeof(Ke));
1018:       for (q=0; q<4; q++) {
1019:         PetscReal   phi[4],dphi[4][2],jw,eta,deta,beta2,dbeta2;
1020:         PetscScalar u,v,du[2],dv[2],h = 0,rbeta2 = 0;
1021:         for (l=0; l<4; l++) {
1022:           phi[l]     = QuadQInterp[q][l];
1023:           dphi[l][0] = QuadQDeriv[q][l][0]*2./hx;
1024:           dphi[l][1] = QuadQDeriv[q][l][1]*2./hy;
1025:           h         += phi[l] * pn[l].h;
1026:           rbeta2    += phi[l] * pn[l].beta2;
1027:         }
1028:         jw = 0.25*hx*hy / thi->rhog; /* rhog is only scaling */
1029:         PointwiseNonlinearity2D(thi,n,phi,dphi,&u,&v,du,dv,&eta,&deta);
1030:         THIFriction(thi,PetscRealPart(rbeta2),PetscRealPart(u*u+v*v)/2,&beta2,&dbeta2);
1031:         for (l=0; l<4; l++) {
1032:           const PetscReal pp = phi[l],*dp = dphi[l];
1033:           for (ll=0; ll<4; ll++) {
1034:             const PetscReal ppl = phi[ll],*dpl = dphi[ll];
1035:             PetscScalar     dgdu,dgdv;
1036:             dgdu = 2.*du[0]*dpl[0] + dv[1]*dpl[0] + 0.5*(du[1]+dv[0])*dpl[1];
1037:             dgdv = 2.*dv[1]*dpl[1] + du[0]*dpl[1] + 0.5*(du[1]+dv[0])*dpl[0];
1038:             /* Picard part */
1039:             Ke[l*2+0][ll*2+0] += dp[0]*jw*eta*4.*dpl[0] + dp[1]*jw*eta*dpl[1] + pp*jw*(beta2/h)*ppl*thi->ssa_friction_scale;
1040:             Ke[l*2+0][ll*2+1] += dp[0]*jw*eta*2.*dpl[1] + dp[1]*jw*eta*dpl[0];
1041:             Ke[l*2+1][ll*2+0] += dp[1]*jw*eta*2.*dpl[0] + dp[0]*jw*eta*dpl[1];
1042:             Ke[l*2+1][ll*2+1] += dp[1]*jw*eta*4.*dpl[1] + dp[0]*jw*eta*dpl[0] + pp*jw*(beta2/h)*ppl*thi->ssa_friction_scale;
1043:             /* extra Newton terms */
1044:             Ke[l*2+0][ll*2+0] += dp[0]*jw*deta*dgdu*(4.*du[0]+2.*dv[1]) + dp[1]*jw*deta*dgdu*(du[1]+dv[0]) + pp*jw*(dbeta2/h)*u*u*ppl*thi->ssa_friction_scale;
1045:             Ke[l*2+0][ll*2+1] += dp[0]*jw*deta*dgdv*(4.*du[0]+2.*dv[1]) + dp[1]*jw*deta*dgdv*(du[1]+dv[0]) + pp*jw*(dbeta2/h)*u*v*ppl*thi->ssa_friction_scale;
1046:             Ke[l*2+1][ll*2+0] += dp[1]*jw*deta*dgdu*(4.*dv[1]+2.*du[0]) + dp[0]*jw*deta*dgdu*(du[1]+dv[0]) + pp*jw*(dbeta2/h)*v*u*ppl*thi->ssa_friction_scale;
1047:             Ke[l*2+1][ll*2+1] += dp[1]*jw*deta*dgdv*(4.*dv[1]+2.*du[0]) + dp[0]*jw*deta*dgdv*(du[1]+dv[0]) + pp*jw*(dbeta2/h)*v*v*ppl*thi->ssa_friction_scale;
1048:           }
1049:         }
1050:       }
1051:       {
1052:         const MatStencil rc[4] = {{0,i,j,0},{0,i+1,j,0},{0,i+1,j+1,0},{0,i,j+1,0}};
1053:         MatSetValuesBlockedStencil(B,4,rc,4,rc,&Ke[0][0],ADD_VALUES);
1054:       }
1055:     }
1056:   }
1057:   THIDARestorePrm(info->da,&prm);

1059:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
1060:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
1061:   MatSetOption(B,MAT_SYMMETRIC,PETSC_TRUE);
1062:   if (thi->verbose) {THIMatrixStatistics(thi,B,PETSC_VIEWER_STDOUT_WORLD);}
1063:   return(0);
1064: }

1068: static PetscErrorCode THIJacobianLocal_3D(DMDALocalInfo *info,Node ***x,Mat B,THI thi,THIAssemblyMode amode)
1069: {
1070:   PetscInt       xs,ys,xm,ym,zm,i,j,k,q,l,ll;
1071:   PetscReal      hx,hy;
1072:   PrmNode        **prm;

1076:   xs = info->zs;
1077:   ys = info->ys;
1078:   xm = info->zm;
1079:   ym = info->ym;
1080:   zm = info->xm;
1081:   hx = thi->Lx / info->mz;
1082:   hy = thi->Ly / info->my;

1084:   MatZeroEntries(B);
1085:   THIDAGetPrm(info->da,&prm);

1087:   for (i=xs; i<xs+xm; i++) {
1088:     for (j=ys; j<ys+ym; j++) {
1089:       PrmNode pn[4];
1090:       QuadExtract(prm,i,j,pn);
1091:       for (k=0; k<zm-1; k++) {
1092:         Node        n[8];
1093:         PetscReal   zn[8],etabase = 0;
1094:         PetscScalar Ke[8*2][8*2];
1095:         PetscInt    ls = 0;

1097:         PrmHexGetZ(pn,k,zm,zn);
1098:         HexExtract(x,i,j,k,n);
1099:         PetscMemzero(Ke,sizeof(Ke));
1100:         if (thi->no_slip && k == 0) {
1101:           for (l=0; l<4; l++) n[l].u = n[l].v = 0;
1102:           ls = 4;
1103:         }
1104:         for (q=0; q<8; q++) {
1105:           PetscReal   dz[3],phi[8],dphi[8][3],jw,eta,deta;
1106:           PetscScalar du[3],dv[3],u,v;
1107:           HexGrad(HexQDeriv[q],zn,dz);
1108:           HexComputeGeometry(q,hx,hy,dz,phi,dphi,&jw);
1109:           PointwiseNonlinearity(thi,n,phi,dphi,&u,&v,du,dv,&eta,&deta);
1110:           jw /= thi->rhog;      /* residuals are scaled by this factor */
1111:           if (q == 0) etabase = eta;
1112:           for (l=ls; l<8; l++) { /* test functions */
1113:             const PetscReal *restrict dp = dphi[l];
1114: #if USE_SSE2_KERNELS
1115:             /* gcc (up to my 4.5 snapshot) is really bad at hoisting intrinsics so we do it manually */
1116:             __m128d
1117:               p4         = _mm_set1_pd(4),p2 = _mm_set1_pd(2),p05 = _mm_set1_pd(0.5),
1118:               p42        = _mm_setr_pd(4,2),p24 = _mm_shuffle_pd(p42,p42,_MM_SHUFFLE2(0,1)),
1119:               du0        = _mm_set1_pd(du[0]),du1 = _mm_set1_pd(du[1]),du2 = _mm_set1_pd(du[2]),
1120:               dv0        = _mm_set1_pd(dv[0]),dv1 = _mm_set1_pd(dv[1]),dv2 = _mm_set1_pd(dv[2]),
1121:               jweta      = _mm_set1_pd(jw*eta),jwdeta = _mm_set1_pd(jw*deta),
1122:               dp0        = _mm_set1_pd(dp[0]),dp1 = _mm_set1_pd(dp[1]),dp2 = _mm_set1_pd(dp[2]),
1123:               dp0jweta   = _mm_mul_pd(dp0,jweta),dp1jweta = _mm_mul_pd(dp1,jweta),dp2jweta = _mm_mul_pd(dp2,jweta),
1124:               p4du0p2dv1 = _mm_add_pd(_mm_mul_pd(p4,du0),_mm_mul_pd(p2,dv1)), /* 4 du0 + 2 dv1 */
1125:               p4dv1p2du0 = _mm_add_pd(_mm_mul_pd(p4,dv1),_mm_mul_pd(p2,du0)), /* 4 dv1 + 2 du0 */
1126:               pdu2dv2    = _mm_unpacklo_pd(du2,dv2),                          /* [du2, dv2] */
1127:               du1pdv0    = _mm_add_pd(du1,dv0),                               /* du1 + dv0 */
1128:               t1         = _mm_mul_pd(dp0,p4du0p2dv1),                        /* dp0 (4 du0 + 2 dv1) */
1129:               t2         = _mm_mul_pd(dp1,p4dv1p2du0);                        /* dp1 (4 dv1 + 2 du0) */

1131: #endif
1132: #if defined COMPUTE_LOWER_TRIANGULAR  /* The element matrices are always symmetric so computing the lower-triangular part is not necessary */
1133:             for (ll=ls; ll<8; ll++) { /* trial functions */
1134: #else
1135:             for (ll=l; ll<8; ll++) {
1136: #endif
1137:               const PetscReal *restrict dpl = dphi[ll];
1138:               if (amode == THIASSEMBLY_TRIDIAGONAL && (l-ll)%4) continue; /* these entries would not be inserted */
1139: #if !USE_SSE2_KERNELS
1140:               /* The analytic Jacobian in nice, easy-to-read form */
1141:               {
1142:                 PetscScalar dgdu,dgdv;
1143:                 dgdu = 2.*du[0]*dpl[0] + dv[1]*dpl[0] + 0.5*(du[1]+dv[0])*dpl[1] + 0.5*du[2]*dpl[2];
1144:                 dgdv = 2.*dv[1]*dpl[1] + du[0]*dpl[1] + 0.5*(du[1]+dv[0])*dpl[0] + 0.5*dv[2]*dpl[2];
1145:                 /* Picard part */
1146:                 Ke[l*2+0][ll*2+0] += dp[0]*jw*eta*4.*dpl[0] + dp[1]*jw*eta*dpl[1] + dp[2]*jw*eta*dpl[2];
1147:                 Ke[l*2+0][ll*2+1] += dp[0]*jw*eta*2.*dpl[1] + dp[1]*jw*eta*dpl[0];
1148:                 Ke[l*2+1][ll*2+0] += dp[1]*jw*eta*2.*dpl[0] + dp[0]*jw*eta*dpl[1];
1149:                 Ke[l*2+1][ll*2+1] += dp[1]*jw*eta*4.*dpl[1] + dp[0]*jw*eta*dpl[0] + dp[2]*jw*eta*dpl[2];
1150:                 /* extra Newton terms */
1151:                 Ke[l*2+0][ll*2+0] += dp[0]*jw*deta*dgdu*(4.*du[0]+2.*dv[1]) + dp[1]*jw*deta*dgdu*(du[1]+dv[0]) + dp[2]*jw*deta*dgdu*du[2];
1152:                 Ke[l*2+0][ll*2+1] += dp[0]*jw*deta*dgdv*(4.*du[0]+2.*dv[1]) + dp[1]*jw*deta*dgdv*(du[1]+dv[0]) + dp[2]*jw*deta*dgdv*du[2];
1153:                 Ke[l*2+1][ll*2+0] += dp[1]*jw*deta*dgdu*(4.*dv[1]+2.*du[0]) + dp[0]*jw*deta*dgdu*(du[1]+dv[0]) + dp[2]*jw*deta*dgdu*dv[2];
1154:                 Ke[l*2+1][ll*2+1] += dp[1]*jw*deta*dgdv*(4.*dv[1]+2.*du[0]) + dp[0]*jw*deta*dgdv*(du[1]+dv[0]) + dp[2]*jw*deta*dgdv*dv[2];
1155:               }
1156: #else
1157:               /* This SSE2 code is an exact replica of above, but uses explicit packed instructions for some speed
1158:               * benefit.  On my hardware, these intrinsics are almost twice as fast as above, reducing total assembly cost
1159:               * by 25 to 30 percent. */
1160:               {
1161:                 __m128d
1162:                   keu   = _mm_loadu_pd(&Ke[l*2+0][ll*2+0]),
1163:                   kev   = _mm_loadu_pd(&Ke[l*2+1][ll*2+0]),
1164:                   dpl01 = _mm_loadu_pd(&dpl[0]),dpl10 = _mm_shuffle_pd(dpl01,dpl01,_MM_SHUFFLE2(0,1)),dpl2 = _mm_set_sd(dpl[2]),
1165:                   t0,t3,pdgduv;
1166:                 keu = _mm_add_pd(keu,_mm_add_pd(_mm_mul_pd(_mm_mul_pd(dp0jweta,p42),dpl01),
1167:                                                 _mm_add_pd(_mm_mul_pd(dp1jweta,dpl10),
1168:                                                            _mm_mul_pd(dp2jweta,dpl2))));
1169:                 kev = _mm_add_pd(kev,_mm_add_pd(_mm_mul_pd(_mm_mul_pd(dp1jweta,p24),dpl01),
1170:                                                 _mm_add_pd(_mm_mul_pd(dp0jweta,dpl10),
1171:                                                            _mm_mul_pd(dp2jweta,_mm_shuffle_pd(dpl2,dpl2,_MM_SHUFFLE2(0,1))))));
1172:                 pdgduv = _mm_mul_pd(p05,_mm_add_pd(_mm_add_pd(_mm_mul_pd(p42,_mm_mul_pd(du0,dpl01)),
1173:                                                               _mm_mul_pd(p24,_mm_mul_pd(dv1,dpl01))),
1174:                                                    _mm_add_pd(_mm_mul_pd(du1pdv0,dpl10),
1175:                                                               _mm_mul_pd(pdu2dv2,_mm_set1_pd(dpl[2]))))); /* [dgdu, dgdv] */
1176:                 t0 = _mm_mul_pd(jwdeta,pdgduv);  /* jw deta [dgdu, dgdv] */
1177:                 t3 = _mm_mul_pd(t0,du1pdv0);     /* t0 (du1 + dv0) */
1178:                 _mm_storeu_pd(&Ke[l*2+0][ll*2+0],_mm_add_pd(keu,_mm_add_pd(_mm_mul_pd(t1,t0),
1179:                                                                            _mm_add_pd(_mm_mul_pd(dp1,t3),
1180:                                                                                       _mm_mul_pd(t0,_mm_mul_pd(dp2,du2))))));
1181:                 _mm_storeu_pd(&Ke[l*2+1][ll*2+0],_mm_add_pd(kev,_mm_add_pd(_mm_mul_pd(t2,t0),
1182:                                                                            _mm_add_pd(_mm_mul_pd(dp0,t3),
1183:                                                                                       _mm_mul_pd(t0,_mm_mul_pd(dp2,dv2))))));
1184:               }
1185: #endif
1186:             }
1187:           }
1188:         }
1189:         if (k == 0) { /* on a bottom face */
1190:           if (thi->no_slip) {
1191:             const PetscReal   hz    = PetscRealPart(pn[0].h)/(zm-1);
1192:             const PetscScalar diagu = 2*etabase/thi->rhog*(hx*hy/hz + hx*hz/hy + 4*hy*hz/hx),diagv = 2*etabase/thi->rhog*(hx*hy/hz + 4*hx*hz/hy + hy*hz/hx);
1193:             Ke[0][0] = thi->dirichlet_scale*diagu;
1194:             Ke[1][1] = thi->dirichlet_scale*diagv;
1195:           } else {
1196:             for (q=0; q<4; q++) {
1197:               const PetscReal jw = 0.25*hx*hy/thi->rhog,*phi = QuadQInterp[q];
1198:               PetscScalar     u  =0,v=0,rbeta2=0;
1199:               PetscReal       beta2,dbeta2;
1200:               for (l=0; l<4; l++) {
1201:                 u      += phi[l]*n[l].u;
1202:                 v      += phi[l]*n[l].v;
1203:                 rbeta2 += phi[l]*pn[l].beta2;
1204:               }
1205:               THIFriction(thi,PetscRealPart(rbeta2),PetscRealPart(u*u+v*v)/2,&beta2,&dbeta2);
1206:               for (l=0; l<4; l++) {
1207:                 const PetscReal pp = phi[l];
1208:                 for (ll=0; ll<4; ll++) {
1209:                   const PetscReal ppl = phi[ll];
1210:                   Ke[l*2+0][ll*2+0] += pp*jw*beta2*ppl + pp*jw*dbeta2*u*u*ppl;
1211:                   Ke[l*2+0][ll*2+1] +=                   pp*jw*dbeta2*u*v*ppl;
1212:                   Ke[l*2+1][ll*2+0] +=                   pp*jw*dbeta2*v*u*ppl;
1213:                   Ke[l*2+1][ll*2+1] += pp*jw*beta2*ppl + pp*jw*dbeta2*v*v*ppl;
1214:                 }
1215:               }
1216:             }
1217:           }
1218:         }
1219:         {
1220:           const MatStencil rc[8] = {{i,j,k,0},{i+1,j,k,0},{i+1,j+1,k,0},{i,j+1,k,0},{i,j,k+1,0},{i+1,j,k+1,0},{i+1,j+1,k+1,0},{i,j+1,k+1,0}};
1221:           if (amode == THIASSEMBLY_TRIDIAGONAL) {
1222:             for (l=0; l<4; l++) { /* Copy out each of the blocks, discarding horizontal coupling */
1223:               const PetscInt   l4     = l+4;
1224:               const MatStencil rcl[2] = {{rc[l].k,rc[l].j,rc[l].i,0},{rc[l4].k,rc[l4].j,rc[l4].i,0}};
1225: #if defined COMPUTE_LOWER_TRIANGULAR
1226:               const PetscScalar Kel[4][4] = {{Ke[2*l+0][2*l+0] ,Ke[2*l+0][2*l+1] ,Ke[2*l+0][2*l4+0] ,Ke[2*l+0][2*l4+1]},
1227:                                              {Ke[2*l+1][2*l+0] ,Ke[2*l+1][2*l+1] ,Ke[2*l+1][2*l4+0] ,Ke[2*l+1][2*l4+1]},
1228:                                              {Ke[2*l4+0][2*l+0],Ke[2*l4+0][2*l+1],Ke[2*l4+0][2*l4+0],Ke[2*l4+0][2*l4+1]},
1229:                                              {Ke[2*l4+1][2*l+0],Ke[2*l4+1][2*l+1],Ke[2*l4+1][2*l4+0],Ke[2*l4+1][2*l4+1]}};
1230: #else
1231:               /* Same as above except for the lower-left block */
1232:               const PetscScalar Kel[4][4] = {{Ke[2*l+0][2*l+0] ,Ke[2*l+0][2*l+1] ,Ke[2*l+0][2*l4+0] ,Ke[2*l+0][2*l4+1]},
1233:                                              {Ke[2*l+1][2*l+0] ,Ke[2*l+1][2*l+1] ,Ke[2*l+1][2*l4+0] ,Ke[2*l+1][2*l4+1]},
1234:                                              {Ke[2*l+0][2*l4+0],Ke[2*l+1][2*l4+0],Ke[2*l4+0][2*l4+0],Ke[2*l4+0][2*l4+1]},
1235:                                              {Ke[2*l+0][2*l4+1],Ke[2*l+1][2*l4+1],Ke[2*l4+1][2*l4+0],Ke[2*l4+1][2*l4+1]}};
1236: #endif
1237:               MatSetValuesBlockedStencil(B,2,rcl,2,rcl,&Kel[0][0],ADD_VALUES);
1238:             }
1239:           } else {
1240: #if !defined COMPUTE_LOWER_TRIANGULAR /* fill in lower-triangular part, this is really cheap compared to computing the entries */
1241:             for (l=0; l<8; l++) {
1242:               for (ll=l+1; ll<8; ll++) {
1243:                 Ke[ll*2+0][l*2+0] = Ke[l*2+0][ll*2+0];
1244:                 Ke[ll*2+1][l*2+0] = Ke[l*2+0][ll*2+1];
1245:                 Ke[ll*2+0][l*2+1] = Ke[l*2+1][ll*2+0];
1246:                 Ke[ll*2+1][l*2+1] = Ke[l*2+1][ll*2+1];
1247:               }
1248:             }
1249: #endif
1250:             MatSetValuesBlockedStencil(B,8,rc,8,rc,&Ke[0][0],ADD_VALUES);
1251:           }
1252:         }
1253:       }
1254:     }
1255:   }
1256:   THIDARestorePrm(info->da,&prm);

1258:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
1259:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
1260:   MatSetOption(B,MAT_SYMMETRIC,PETSC_TRUE);
1261:   if (thi->verbose) {THIMatrixStatistics(thi,B,PETSC_VIEWER_STDOUT_WORLD);}
1262:   return(0);
1263: }

1267: static PetscErrorCode THIJacobianLocal_3D_Full(DMDALocalInfo *info,Node ***x,Mat A,Mat B,THI thi)
1268: {

1272:   THIJacobianLocal_3D(info,x,B,thi,THIASSEMBLY_FULL);
1273:   return(0);
1274: }

1278: static PetscErrorCode THIJacobianLocal_3D_Tridiagonal(DMDALocalInfo *info,Node ***x,Mat B,THI thi)
1279: {

1283:   THIJacobianLocal_3D(info,x,B,thi,THIASSEMBLY_TRIDIAGONAL);
1284:   return(0);
1285: }

1289: static PetscErrorCode DMRefineHierarchy_THI(DM dac0,PetscInt nlevels,DM hierarchy[])
1290: {
1291:   PetscErrorCode  ierr;
1292:   THI             thi;
1293:   PetscInt        dim,M,N,m,n,s,dof;
1294:   DM              dac,daf;
1295:   DMDAStencilType st;
1296:   DM_DA           *ddf,*ddc;

1299:   PetscObjectQuery((PetscObject)dac0,"THI",(PetscObject*)&thi);
1300:   if (!thi) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Cannot refine this DMDA, missing composed THI instance");
1301:   if (nlevels > 1) {
1302:     DMRefineHierarchy(dac0,nlevels-1,hierarchy);
1303:     dac  = hierarchy[nlevels-2];
1304:   } else {
1305:     dac = dac0;
1306:   }
1307:   DMDAGetInfo(dac,&dim, &N,&M,0, &n,&m,0, &dof,&s,0,0,0,&st);
1308:   if (dim != 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"This function can only refine 2D DMDAs");

1310:   /* Creates a 3D DMDA with the same map-plane layout as the 2D one, with contiguous columns */
1311:   DMDACreate3d(PetscObjectComm((PetscObject)dac),DM_BOUNDARY_NONE,DM_BOUNDARY_PERIODIC,DM_BOUNDARY_PERIODIC,st,thi->zlevels,N,M,1,n,m,dof,s,NULL,NULL,NULL,&daf);

1313:   daf->ops->creatematrix        = dac->ops->creatematrix;
1314:   daf->ops->createinterpolation = dac->ops->createinterpolation;
1315:   daf->ops->getcoloring         = dac->ops->getcoloring;
1316:   ddf                           = (DM_DA*)daf->data;
1317:   ddc                           = (DM_DA*)dac->data;
1318:   ddf->interptype               = ddc->interptype;

1320:   DMDASetFieldName(daf,0,"x-velocity");
1321:   DMDASetFieldName(daf,1,"y-velocity");

1323:   hierarchy[nlevels-1] = daf;
1324:   return(0);
1325: }

1329: static PetscErrorCode DMCreateInterpolation_DA_THI(DM dac,DM daf,Mat *A,Vec *scale)
1330: {
1332:   PetscInt       dim;

1339:   DMDAGetInfo(daf,&dim,0,0,0,0,0,0,0,0,0,0,0,0);
1340:   if (dim  == 2) {
1341:     /* We are in the 2D problem and use normal DMDA interpolation */
1342:     DMCreateInterpolation(dac,daf,A,scale);
1343:   } else {
1344:     PetscInt i,j,k,xs,ys,zs,xm,ym,zm,mx,my,mz,rstart,cstart;
1345:     Mat      B;

1347:     DMDAGetInfo(daf,0, &mz,&my,&mx, 0,0,0, 0,0,0,0,0,0);
1348:     DMDAGetCorners(daf,&zs,&ys,&xs,&zm,&ym,&xm);
1349:     if (zs != 0) SETERRQ(PETSC_COMM_SELF,1,"unexpected");
1350:     MatCreate(PetscObjectComm((PetscObject)daf),&B);
1351:     MatSetSizes(B,xm*ym*zm,xm*ym,mx*my*mz,mx*my);

1353:     MatSetType(B,MATAIJ);
1354:     MatSeqAIJSetPreallocation(B,1,NULL);
1355:     MatMPIAIJSetPreallocation(B,1,NULL,0,NULL);
1356:     MatGetOwnershipRange(B,&rstart,NULL);
1357:     MatGetOwnershipRangeColumn(B,&cstart,NULL);
1358:     for (i=xs; i<xs+xm; i++) {
1359:       for (j=ys; j<ys+ym; j++) {
1360:         for (k=zs; k<zs+zm; k++) {
1361:           PetscInt    i2  = i*ym+j,i3 = i2*zm+k;
1362:           PetscScalar val = ((k == 0 || k == mz-1) ? 0.5 : 1.) / (mz-1.); /* Integration using trapezoid rule */
1363:           MatSetValue(B,cstart+i3,rstart+i2,val,INSERT_VALUES);
1364:         }
1365:       }
1366:     }
1367:     MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
1368:     MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
1369:     MatCreateMAIJ(B,sizeof(Node)/sizeof(PetscScalar),A);
1370:     MatDestroy(&B);
1371:   }
1372:   return(0);
1373: }

1377: static PetscErrorCode DMCreateMatrix_THI_Tridiagonal(DM da,Mat *J)
1378: {
1379:   PetscErrorCode         ierr;
1380:   Mat                    A;
1381:   PetscInt               xm,ym,zm,dim,dof = 2,starts[3],dims[3];
1382:   ISLocalToGlobalMapping ltog,ltogb;

1385:   DMDAGetInfo(da,&dim, 0,0,0, 0,0,0, 0,0,0,0,0,0);
1386:   if (dim != 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected DMDA to be 3D");
1387:   DMDAGetCorners(da,0,0,0,&zm,&ym,&xm);
1388:   DMGetLocalToGlobalMapping(da,&ltog);
1389:   DMGetLocalToGlobalMappingBlock(da,&ltogb);
1390:   MatCreate(PetscObjectComm((PetscObject)da),&A);
1391:   MatSetSizes(A,dof*xm*ym*zm,dof*xm*ym*zm,PETSC_DETERMINE,PETSC_DETERMINE);
1392:   MatSetType(A,da->mattype);
1393:   MatSetFromOptions(A);
1394:   MatSeqAIJSetPreallocation(A,3*2,NULL);
1395:   MatMPIAIJSetPreallocation(A,3*2,NULL,0,NULL);
1396:   MatSeqBAIJSetPreallocation(A,2,3,NULL);
1397:   MatMPIBAIJSetPreallocation(A,2,3,NULL,0,NULL);
1398:   MatSeqSBAIJSetPreallocation(A,2,2,NULL);
1399:   MatMPISBAIJSetPreallocation(A,2,2,NULL,0,NULL);
1400:   MatSetLocalToGlobalMapping(A,ltog,ltog);
1401:   MatSetLocalToGlobalMappingBlock(A,ltogb,ltogb);
1402:   DMDAGetGhostCorners(da,&starts[0],&starts[1],&starts[2],&dims[0],&dims[1],&dims[2]);
1403:   MatSetStencil(A,dim,dims,starts,dof);
1404:   *J   = A;
1405:   return(0);
1406: }

1410: static PetscErrorCode THIDAVecView_VTK_XML(THI thi,DM da,Vec X,const char filename[])
1411: {
1412:   const PetscInt dof   = 2;
1413:   Units          units = thi->units;
1414:   MPI_Comm       comm;
1416:   PetscViewer    viewer;
1417:   PetscMPIInt    rank,size,tag,nn,nmax;
1418:   PetscInt       mx,my,mz,r,range[6];
1419:   PetscScalar    *x;

1422:   PetscObjectGetComm((PetscObject)thi,&comm);
1423:   DMDAGetInfo(da,0, &mz,&my,&mx, 0,0,0, 0,0,0,0,0,0);
1424:   MPI_Comm_size(comm,&size);
1425:   MPI_Comm_rank(comm,&rank);
1426:   PetscViewerASCIIOpen(comm,filename,&viewer);
1427:   PetscViewerASCIIPrintf(viewer,"<VTKFile type=\"StructuredGrid\" version=\"0.1\" byte_order=\"LittleEndian\">\n");
1428:   PetscViewerASCIIPrintf(viewer,"  <StructuredGrid WholeExtent=\"%d %D %d %D %d %D\">\n",0,mz-1,0,my-1,0,mx-1);

1430:   DMDAGetCorners(da,range,range+1,range+2,range+3,range+4,range+5);
1431:   PetscMPIIntCast(range[3]*range[4]*range[5]*dof,&nn);
1432:   MPI_Reduce(&nn,&nmax,1,MPI_INT,MPI_MAX,0,comm);
1433:   tag  = ((PetscObject) viewer)->tag;
1434:   VecGetArray(X,&x);
1435:   if (!rank) {
1436:     PetscScalar *array;
1437:     PetscMalloc1(nmax,&array);
1438:     for (r=0; r<size; r++) {
1439:       PetscInt    i,j,k,xs,xm,ys,ym,zs,zm;
1440:       PetscScalar *ptr;
1441:       MPI_Status  status;
1442:       if (r) {
1443:         MPI_Recv(range,6,MPIU_INT,r,tag,comm,MPI_STATUS_IGNORE);
1444:       }
1445:       zs = range[0];ys = range[1];xs = range[2];zm = range[3];ym = range[4];xm = range[5];
1446:       if (xm*ym*zm*dof > nmax) SETERRQ(PETSC_COMM_SELF,1,"should not happen");
1447:       if (r) {
1448:         MPI_Recv(array,nmax,MPIU_SCALAR,r,tag,comm,&status);
1449:         MPI_Get_count(&status,MPIU_SCALAR,&nn);
1450:         if (nn != xm*ym*zm*dof) SETERRQ(PETSC_COMM_SELF,1,"should not happen");
1451:         ptr = array;
1452:       } else ptr = x;
1453:       PetscViewerASCIIPrintf(viewer,"    <Piece Extent=\"%D %D %D %D %D %D\">\n",zs,zs+zm-1,ys,ys+ym-1,xs,xs+xm-1);

1455:       PetscViewerASCIIPrintf(viewer,"      <Points>\n");
1456:       PetscViewerASCIIPrintf(viewer,"        <DataArray type=\"Float32\" NumberOfComponents=\"3\" format=\"ascii\">\n");
1457:       for (i=xs; i<xs+xm; i++) {
1458:         for (j=ys; j<ys+ym; j++) {
1459:           for (k=zs; k<zs+zm; k++) {
1460:             PrmNode   p;
1461:             PetscReal xx = thi->Lx*i/mx,yy = thi->Ly*j/my,zz;
1462:             thi->initialize(thi,xx,yy,&p);
1463:             zz   = PetscRealPart(p.b) + PetscRealPart(p.h)*k/(mz-1);
1464:             PetscViewerASCIIPrintf(viewer,"%f %f %f\n",(double)xx,(double)yy,(double)zz);
1465:           }
1466:         }
1467:       }
1468:       PetscViewerASCIIPrintf(viewer,"        </DataArray>\n");
1469:       PetscViewerASCIIPrintf(viewer,"      </Points>\n");

1471:       PetscViewerASCIIPrintf(viewer,"      <PointData>\n");
1472:       PetscViewerASCIIPrintf(viewer,"        <DataArray type=\"Float32\" Name=\"velocity\" NumberOfComponents=\"3\" format=\"ascii\">\n");
1473:       for (i=0; i<nn; i+=dof) {
1474:         PetscViewerASCIIPrintf(viewer,"%f %f %f\n",(double)(PetscRealPart(ptr[i])*units->year/units->meter),(double)(PetscRealPart(ptr[i+1])*units->year/units->meter),0.0);
1475:       }
1476:       PetscViewerASCIIPrintf(viewer,"        </DataArray>\n");

1478:       PetscViewerASCIIPrintf(viewer,"        <DataArray type=\"Int32\" Name=\"rank\" NumberOfComponents=\"1\" format=\"ascii\">\n");
1479:       for (i=0; i<nn; i+=dof) {
1480:         PetscViewerASCIIPrintf(viewer,"%D\n",r);
1481:       }
1482:       PetscViewerASCIIPrintf(viewer,"        </DataArray>\n");
1483:       PetscViewerASCIIPrintf(viewer,"      </PointData>\n");

1485:       PetscViewerASCIIPrintf(viewer,"    </Piece>\n");
1486:     }
1487:     PetscFree(array);
1488:   } else {
1489:     MPI_Send(range,6,MPIU_INT,0,tag,comm);
1490:     MPI_Send(x,nn,MPIU_SCALAR,0,tag,comm);
1491:   }
1492:   VecRestoreArray(X,&x);
1493:   PetscViewerASCIIPrintf(viewer,"  </StructuredGrid>\n");
1494:   PetscViewerASCIIPrintf(viewer,"</VTKFile>\n");
1495:   PetscViewerDestroy(&viewer);
1496:   return(0);
1497: }

1501: int main(int argc,char *argv[])
1502: {
1503:   MPI_Comm       comm;
1504:   THI            thi;
1506:   DM             da;
1507:   SNES           snes;

1509:   PetscInitialize(&argc,&argv,0,help);
1510:   comm = PETSC_COMM_WORLD;

1512:   THICreate(comm,&thi);
1513:   {
1514:     PetscInt M = 3,N = 3,P = 2;
1515:     PetscOptionsBegin(comm,NULL,"Grid resolution options","");
1516:     {
1517:       PetscOptionsInt("-M","Number of elements in x-direction on coarse level","",M,&M,NULL);
1518:       N    = M;
1519:       PetscOptionsInt("-N","Number of elements in y-direction on coarse level (if different from M)","",N,&N,NULL);
1520:       if (thi->coarse2d) {
1521:         PetscOptionsInt("-zlevels","Number of elements in z-direction on fine level","",thi->zlevels,&thi->zlevels,NULL);
1522:       } else {
1523:         PetscOptionsInt("-P","Number of elements in z-direction on coarse level","",P,&P,NULL);
1524:       }
1525:     }
1526:     PetscOptionsEnd();
1527:     if (thi->coarse2d) {
1528:       DMDACreate2d(comm,DM_BOUNDARY_PERIODIC,DM_BOUNDARY_PERIODIC,DMDA_STENCIL_BOX,-N,-M,PETSC_DETERMINE,PETSC_DETERMINE,sizeof(Node)/sizeof(PetscScalar),1,0,0,&da);

1530:       da->ops->refinehierarchy     = DMRefineHierarchy_THI;
1531:       da->ops->createinterpolation = DMCreateInterpolation_DA_THI;

1533:       PetscObjectCompose((PetscObject)da,"THI",(PetscObject)thi);
1534:     } else {
1535:       DMDACreate3d(comm,DM_BOUNDARY_NONE,DM_BOUNDARY_PERIODIC,DM_BOUNDARY_PERIODIC, DMDA_STENCIL_BOX,-P,-N,-M,1,PETSC_DETERMINE,PETSC_DETERMINE,sizeof(Node)/sizeof(PetscScalar),1,0,0,0,&da);
1536:     }
1537:     DMDASetFieldName(da,0,"x-velocity");
1538:     DMDASetFieldName(da,1,"y-velocity");
1539:   }
1540:   THISetUpDM(thi,da);
1541:   if (thi->tridiagonal) da->ops->creatematrix = DMCreateMatrix_THI_Tridiagonal;

1543:   {                             /* Set the fine level matrix type if -da_refine */
1544:     PetscInt rlevel,clevel;
1545:     DMGetRefineLevel(da,&rlevel);
1546:     DMGetCoarsenLevel(da,&clevel);
1547:     if (rlevel - clevel > 0) {DMSetMatType(da,thi->mattype);}
1548:   }

1550:   DMDASNESSetFunctionLocal(da,ADD_VALUES,(DMDASNESFunction)THIFunctionLocal,thi);
1551:   if (thi->tridiagonal) {
1552:     DMDASNESSetJacobianLocal(da,(DMDASNESJacobian)THIJacobianLocal_3D_Tridiagonal,thi);
1553:   } else {
1554:     DMDASNESSetJacobianLocal(da,(DMDASNESJacobian)THIJacobianLocal_3D_Full,thi);
1555:   }
1556:   DMCoarsenHookAdd(da,DMCoarsenHook_THI,NULL,thi);
1557:   DMRefineHookAdd(da,DMRefineHook_THI,NULL,thi);

1559:   DMSetApplicationContext(da,thi);

1561:   SNESCreate(comm,&snes);
1562:   SNESSetDM(snes,da);
1563:   DMDestroy(&da);
1564:   SNESSetComputeInitialGuess(snes,THIInitial,NULL);
1565:   SNESSetFromOptions(snes);

1567:   SNESSolve(snes,NULL,NULL);

1569:   THISolveStatistics(thi,snes,0,"Full");

1571:   {
1572:     PetscBool flg;
1573:     char      filename[PETSC_MAX_PATH_LEN] = "";
1574:     PetscOptionsGetString(NULL,"-o",filename,sizeof(filename),&flg);
1575:     if (flg) {
1576:       Vec X;
1577:       DM  dm;
1578:       SNESGetSolution(snes,&X);
1579:       SNESGetDM(snes,&dm);
1580:       THIDAVecView_VTK_XML(thi,dm,X,filename);
1581:     }
1582:   }

1584:   DMDestroy(&da);
1585:   SNESDestroy(&snes);
1586:   THIDestroy(&thi);
1587:   PetscFinalize();
1588:   return 0;
1589: }