Actual source code: ex48.c

petsc-3.5.4 2015-05-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: #if defined(PETSC_APPLE_FRAMEWORK)
 59: #import <PETSc/petscsnes.h>
 60: #import <PETSc/petsc-private/dmdaimpl.h>     /* There is not yet a public interface to manipulate dm->ops */
 61: #else
 62: #include <petscsnes.h>
 63: #include <petsc-private/dmdaimpl.h>     /* There is not yet a public interface to manipulate dm->ops */
 64: #endif
 65: #include <ctype.h>              /* toupper() */

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

 78: /* The SSE2 kernels are only for PetscScalar=double on architectures that support it */
 79: #define USE_SSE2_KERNELS (!defined NO_SSE2                              \
 80:                           && !defined PETSC_USE_COMPLEX                 \
 81:                           && !defined PETSC_USE_REAL_SINGLE             \
 82:                           && !defined PETSC_USE_REAL___FLOAT128         \
 83:                           && defined __SSE2__)

 85: static PetscClassId THI_CLASSID;

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

149: #define HexExtract(x,i,j,k,n) do {              \
150:     (n)[0] = (x)[i][j][k];                      \
151:     (n)[1] = (x)[i+1][j][k];                    \
152:     (n)[2] = (x)[i+1][j+1][k];                  \
153:     (n)[3] = (x)[i][j+1][k];                    \
154:     (n)[4] = (x)[i][j][k+1];                    \
155:     (n)[5] = (x)[i+1][j][k+1];                  \
156:     (n)[6] = (x)[i+1][j+1][k+1];                \
157:     (n)[7] = (x)[i][j+1][k+1];                  \
158:   } while (0)

160: #define HexExtractRef(x,i,j,k,n) do {           \
161:     (n)[0] = &(x)[i][j][k];                     \
162:     (n)[1] = &(x)[i+1][j][k];                   \
163:     (n)[2] = &(x)[i+1][j+1][k];                 \
164:     (n)[3] = &(x)[i][j+1][k];                   \
165:     (n)[4] = &(x)[i][j][k+1];                   \
166:     (n)[5] = &(x)[i+1][j][k+1];                 \
167:     (n)[6] = &(x)[i+1][j+1][k+1];               \
168:     (n)[7] = &(x)[i][j+1][k+1];                 \
169:   } while (0)

171: #define QuadExtract(x,i,j,n) do {               \
172:     (n)[0] = (x)[i][j];                         \
173:     (n)[1] = (x)[i+1][j];                       \
174:     (n)[2] = (x)[i+1][j+1];                     \
175:     (n)[3] = (x)[i][j+1];                       \
176:   } while (0)

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

180: static void HexGrad(const PetscReal dphi[][3],const PetscReal zn[],PetscReal dz[])
181: {
182:   PetscInt i;
183:   dz[0] = dz[1] = dz[2] = 0;
184:   for (i=0; i<8; i++) {
185:     dz[0] += dphi[i][0] * zn[i];
186:     dz[1] += dphi[i][1] * zn[i];
187:     dz[2] += dphi[i][2] * zn[i];
188:   }
189: }

191: static void HexComputeGeometry(PetscInt q,PetscReal hx,PetscReal hy,const PetscReal dz[restrict],PetscReal phi[restrict],PetscReal dphi[restrict][3],PetscReal *restrict jw)
192: {
193:   const PetscReal jac[3][3]  = {{hx/2,0,0}, {0,hy/2,0}, {dz[0],dz[1],dz[2]}};
194:   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]}};
195:   const PetscReal jdet       = jac[0][0]*jac[1][1]*jac[2][2];
196:   PetscInt        i;

198:   for (i=0; i<8; i++) {
199:     const PetscReal *dphir = HexQDeriv[q][i];
200:     phi[i]     = HexQInterp[q][i];
201:     dphi[i][0] = dphir[0]*ijac[0][0] + dphir[1]*ijac[1][0] + dphir[2]*ijac[2][0];
202:     dphi[i][1] = dphir[0]*ijac[0][1] + dphir[1]*ijac[1][1] + dphir[2]*ijac[2][1];
203:     dphi[i][2] = dphir[0]*ijac[0][2] + dphir[1]*ijac[1][2] + dphir[2]*ijac[2][2];
204:   }
205:   *jw = 1.0 * jdet;
206: }

208: typedef struct _p_THI   *THI;
209: typedef struct _n_Units *Units;

211: typedef struct {
212:   PetscScalar u,v;
213: } Node;

215: typedef struct {
216:   PetscScalar b;                /* bed */
217:   PetscScalar h;                /* thickness */
218:   PetscScalar beta2;            /* friction */
219: } PrmNode;

221: typedef struct {
222:   PetscReal min,max,cmin,cmax;
223: } PRange;

225: typedef enum {THIASSEMBLY_TRIDIAGONAL,THIASSEMBLY_FULL} THIAssemblyMode;

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

252: struct _n_Units {
253:   /* fundamental */
254:   PetscReal meter;
255:   PetscReal kilogram;
256:   PetscReal second;
257:   /* derived */
258:   PetscReal Pascal;
259:   PetscReal year;
260: };

262: static PetscErrorCode THIJacobianLocal_3D_Full(DMDALocalInfo*,Node***,Mat,Mat,THI);
263: static PetscErrorCode THIJacobianLocal_3D_Tridiagonal(DMDALocalInfo*,Node***,Mat,THI);
264: static PetscErrorCode THIJacobianLocal_2D(DMDALocalInfo*,Node**,Mat,THI);

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

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

287:   p->b     = s - 1000*units->meter + 500*units->meter * PetscSinReal(x*2*PETSC_PI/thi->Lx) * PetscSinReal(y*2*PETSC_PI/thi->Ly);
288:   p->h     = s - p->b;
289:   p->beta2 = 1e30;
290: }

292: static void THIInitialize_HOM_C(THI thi,PetscReal x,PetscReal y,PrmNode *p)
293: {
294:   Units     units = thi->units;
295:   PetscReal s     = -x*PetscSinReal(thi->alpha);

297:   p->b = s - 1000*units->meter;
298:   p->h = s - p->b;
299:   /* tau_b = beta2 v   is a stress (Pa) */
300:   p->beta2 = 1000 * (1 + PetscSinReal(x*2*PETSC_PI/thi->Lx)*PetscSinReal(y*2*PETSC_PI/thi->Ly)) * units->Pascal * units->year / units->meter;
301: }

303: /* These are just toys */

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

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

323:   p->b = s - 1000*units->meter + 500*units->meter * PetscSinReal(x + PETSC_PI) * PetscSinReal(y + PETSC_PI);
324:   if (PetscRealPart(p->b) > -700*units->meter) p->b += 200*units->meter;
325:   p->h     = s - p->b;
326:   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;
327: }

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

336:   p->b     = s - 1000*units->meter + 500*units->meter * PetscSinReal(x + PETSC_PI) * PetscSinReal(y + PETSC_PI);
337:   p->h     = s - p->b;
338:   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;
339: }

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

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

378: static void RangeUpdate(PetscReal *min,PetscReal *max,PetscReal x)
379: {
380:   if (x < *min) *min = x;
381:   if (x > *max) *max = x;
382: }

384: static void PRangeClear(PRange *p)
385: {
386:   p->cmin = p->min = 1e100;
387:   p->cmax = p->max = -1e100;
388: }

392: static PetscErrorCode PRangeMinMax(PRange *p,PetscReal min,PetscReal max)
393: {

396:   p->cmin = min;
397:   p->cmax = max;
398:   if (min < p->min) p->min = min;
399:   if (max > p->max) p->max = max;
400:   return(0);
401: }

405: static PetscErrorCode THIDestroy(THI *thi)
406: {

410:   if (!*thi) return(0);
411:   if (--((PetscObject)(*thi))->refct > 0) {*thi = 0; return(0);}
412:   PetscFree((*thi)->units);
413:   PetscFree((*thi)->mattype);
414:   PetscHeaderDestroy(thi);
415:   return(0);
416: }

420: static PetscErrorCode THICreate(MPI_Comm comm,THI *inthi)
421: {
422:   static PetscBool registered = PETSC_FALSE;
423:   THI              thi;
424:   Units            units;
425:   PetscErrorCode   ierr;

428:   *inthi = 0;
429:   if (!registered) {
430:     PetscClassIdRegister("Toy Hydrostatic Ice",&THI_CLASSID);
431:     registered = PETSC_TRUE;
432:   }
433:   PetscHeaderCreate(thi,_p_THI,0,THI_CLASSID,"THI","Toy Hydrostatic Ice","",comm,THIDestroy,0);

435:   PetscNew(&thi->units);
436:   units           = thi->units;
437:   units->meter    = 1e-2;
438:   units->second   = 1e-7;
439:   units->kilogram = 1e-12;

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

451:   thi->Lx              = 10.e3;
452:   thi->Ly              = 10.e3;
453:   thi->Lz              = 1000;
454:   thi->dirichlet_scale = 1;
455:   thi->verbose         = PETSC_FALSE;

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

513:     thi->friction.refvel = 100.;
514:     thi->friction.epsvel = 1.;

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

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

522:     PetscOptionsReal("-thi_dirichlet_scale","Scale Dirichlet boundary conditions by this factor","",thi->dirichlet_scale,&thi->dirichlet_scale,NULL);
523:     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);
524:     PetscOptionsBool("-thi_coarse2d","Use a 2D coarse space corresponding to SSA","",thi->coarse2d,&thi->coarse2d,NULL);
525:     PetscOptionsBool("-thi_tridiagonal","Assemble a tridiagonal system (column coupling only) on the finest level","",thi->tridiagonal,&thi->tridiagonal,NULL);
526:     PetscOptionsFList("-thi_mat_type","Matrix type","MatSetType",MatList,mtype,(char*)mtype,sizeof(mtype),NULL);
527:     PetscStrallocpy(mtype,(char**)&thi->mattype);
528:     PetscOptionsBool("-thi_verbose","Enable verbose output (like matrix sizes and statistics)","",thi->verbose,&thi->verbose,NULL);
529:   }
530:   PetscOptionsEnd();

532:   /* dimensionalize */
533:   thi->Lx    *= units->meter;
534:   thi->Ly    *= units->meter;
535:   thi->Lz    *= units->meter;
536:   thi->alpha *= PETSC_PI / 180;

538:   PRangeClear(&thi->eta);
539:   PRangeClear(&thi->beta2);

541:   {
542:     PetscReal u       = 1000*units->meter/(3e7*units->second),
543:               gradu   = u / (100*units->meter),eta,deta,
544:               rho     = 910 * units->kilogram/PetscPowReal(units->meter,3),
545:               grav    = 9.81 * units->meter/PetscSqr(units->second),
546:               driving = rho * grav * PetscSinReal(thi->alpha) * 1000*units->meter;
547:     THIViscosity(thi,0.5*gradu*gradu,&eta,&deta);
548:     thi->rhog = rho * grav;
549:     if (thi->verbose) {
550:       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);
551:       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);
552:       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));
553:       THIViscosity(thi,0.5*PetscSqr(1e-3*gradu),&eta,&deta);
554:       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));
555:     }
556:   }

558:   *inthi = thi;
559:   return(0);
560: }

564: static PetscErrorCode THIInitializePrm(THI thi,DM da2prm,Vec prm)
565: {
566:   PrmNode        **p;
567:   PetscInt       i,j,xs,xm,ys,ym,mx,my;

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

586: static PetscErrorCode THISetUpDM(THI thi,DM dm)
587: {
588:   PetscErrorCode  ierr;
589:   PetscInt        refinelevel,coarsenlevel,level,dim,Mx,My,Mz,mx,my,s;
590:   DMDAStencilType st;
591:   DM              da2prm;
592:   Vec             X;

595:   DMDAGetInfo(dm,&dim, &Mz,&My,&Mx, 0,&my,&mx, 0,&s,0,0,0,&st);
596:   if (dim == 2) {
597:     DMDAGetInfo(dm,&dim, &My,&Mx,0, &my,&mx,0, 0,&s,0,0,0,&st);
598:   }
599:   DMGetRefineLevel(dm,&refinelevel);
600:   DMGetCoarsenLevel(dm,&coarsenlevel);
601:   level = refinelevel - coarsenlevel;
602:   DMDACreate2d(PetscObjectComm((PetscObject)thi),DM_BOUNDARY_PERIODIC,DM_BOUNDARY_PERIODIC,st,My,Mx,my,mx,sizeof(PrmNode)/sizeof(PetscScalar),s,0,0,&da2prm);
603:   DMCreateLocalVector(da2prm,&X);
604:   {
605:     PetscReal Lx = thi->Lx / thi->units->meter,Ly = thi->Ly / thi->units->meter,Lz = thi->Lz / thi->units->meter;
606:     if (dim == 2) {
607:       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));
608:     } else {
609:       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)));
610:     }
611:   }
612:   THIInitializePrm(thi,da2prm,X);
613:   if (thi->tridiagonal) {       /* Reset coarse Jacobian evaluation */
614:     DMDASNESSetJacobianLocal(dm,(DMDASNESJacobian)THIJacobianLocal_3D_Full,thi);
615:   }
616:   if (thi->coarse2d) {
617:     DMDASNESSetJacobianLocal(dm,(DMDASNESJacobian)THIJacobianLocal_2D,thi);
618:   }
619:   PetscObjectCompose((PetscObject)dm,"DMDA2Prm",(PetscObject)da2prm);
620:   PetscObjectCompose((PetscObject)dm,"DMDA2Prm_Vec",(PetscObject)X);
621:   DMDestroy(&da2prm);
622:   VecDestroy(&X);
623:   return(0);
624: }

628: static PetscErrorCode DMCoarsenHook_THI(DM dmf,DM dmc,void *ctx)
629: {
630:   THI            thi = (THI)ctx;
632:   PetscInt       rlevel,clevel;

635:   THISetUpDM(thi,dmc);
636:   DMGetRefineLevel(dmc,&rlevel);
637:   DMGetCoarsenLevel(dmc,&clevel);
638:   if (rlevel-clevel == 0) {DMSetMatType(dmc,MATAIJ);}
639:   DMCoarsenHookAdd(dmc,DMCoarsenHook_THI,NULL,thi);
640:   return(0);
641: }

645: static PetscErrorCode DMRefineHook_THI(DM dmc,DM dmf,void *ctx)
646: {
647:   THI            thi = (THI)ctx;

651:   THISetUpDM(thi,dmf);
652:   DMSetMatType(dmf,thi->mattype);
653:   DMRefineHookAdd(dmf,DMRefineHook_THI,NULL,thi);
654:   /* With grid sequencing, a formerly-refined DM will later be coarsened by PCSetUp_MG */
655:   DMCoarsenHookAdd(dmf,DMCoarsenHook_THI,NULL,thi);
656:   return(0);
657: }

661: static PetscErrorCode THIDAGetPrm(DM da,PrmNode ***prm)
662: {
664:   DM             da2prm;
665:   Vec            X;

668:   PetscObjectQuery((PetscObject)da,"DMDA2Prm",(PetscObject*)&da2prm);
669:   if (!da2prm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"No DMDA2Prm composed with given DMDA");
670:   PetscObjectQuery((PetscObject)da,"DMDA2Prm_Vec",(PetscObject*)&X);
671:   if (!X) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"No DMDA2Prm_Vec composed with given DMDA");
672:   DMDAVecGetArray(da2prm,X,prm);
673:   return(0);
674: }

678: static PetscErrorCode THIDARestorePrm(DM da,PrmNode ***prm)
679: {
681:   DM             da2prm;
682:   Vec            X;

685:   PetscObjectQuery((PetscObject)da,"DMDA2Prm",(PetscObject*)&da2prm);
686:   if (!da2prm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"No DMDA2Prm composed with given DMDA");
687:   PetscObjectQuery((PetscObject)da,"DMDA2Prm_Vec",(PetscObject*)&X);
688:   if (!X) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"No DMDA2Prm_Vec composed with given DMDA");
689:   DMDAVecRestoreArray(da2prm,X,prm);
690:   return(0);
691: }

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

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

730: 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)
731: {
732:   PetscInt    l,ll;
733:   PetscScalar gam;

735:   du[0] = du[1] = du[2] = 0;
736:   dv[0] = dv[1] = dv[2] = 0;
737:   *u    = 0;
738:   *v    = 0;
739:   for (l=0; l<8; l++) {
740:     *u += phi[l] * n[l].u;
741:     *v += phi[l] * n[l].v;
742:     for (ll=0; ll<3; ll++) {
743:       du[ll] += dphi[l][ll] * n[l].u;
744:       dv[ll] += dphi[l][ll] * n[l].v;
745:     }
746:   }
747:   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]);
748:   THIViscosity(thi,PetscRealPart(gam),eta,deta);
749: }

751: 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)
752: {
753:   PetscInt    l,ll;
754:   PetscScalar gam;

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

774: static PetscErrorCode THIFunctionLocal(DMDALocalInfo *info,Node ***x,Node ***f,THI thi)
775: {
776:   PetscInt       xs,ys,xm,ym,zm,i,j,k,q,l;
777:   PetscReal      hx,hy,etamin,etamax,beta2min,beta2max;
778:   PrmNode        **prm;

782:   xs = info->zs;
783:   ys = info->ys;
784:   xm = info->zm;
785:   ym = info->ym;
786:   zm = info->xm;
787:   hx = thi->Lx / info->mz;
788:   hy = thi->Ly / info->my;

790:   etamin   = 1e100;
791:   etamax   = 0;
792:   beta2min = 1e100;
793:   beta2max = 0;

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

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

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

874:   PRangeMinMax(&thi->eta,etamin,etamax);
875:   PRangeMinMax(&thi->beta2,beta2min,beta2max);
876:   return(0);
877: }

881: static PetscErrorCode THIMatrixStatistics(THI thi,Mat B,PetscViewer viewer)
882: {
884:   PetscReal      nrm;
885:   PetscInt       m;
886:   PetscMPIInt    rank;

889:   MatNorm(B,NORM_FROBENIUS,&nrm);
890:   MatGetSize(B,&m,0);
891:   MPI_Comm_rank(PetscObjectComm((PetscObject)B),&rank);
892:   if (!rank) {
893:     PetscScalar val0,val2;
894:     MatGetValue(B,0,0,&val0);
895:     MatGetValue(B,2,2,&val2);
896:     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);
897:   }
898:   return(0);
899: }

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

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

934: static PetscErrorCode THISolveStatistics(THI thi,SNES snes,PetscInt coarsened,const char name[])
935: {
936:   MPI_Comm       comm;
937:   Vec            X;
938:   DM             dm;

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

997: static PetscErrorCode THIJacobianLocal_2D(DMDALocalInfo *info,Node **x,Mat B,THI thi)
998: {
999:   PetscInt       xs,ys,xm,ym,i,j,q,l,ll;
1000:   PetscReal      hx,hy;
1001:   PrmNode        **prm;

1005:   xs = info->ys;
1006:   ys = info->xs;
1007:   xm = info->ym;
1008:   ym = info->xm;
1009:   hx = thi->Lx / info->my;
1010:   hy = thi->Ly / info->mx;

1012:   MatZeroEntries(B);
1013:   THIDAGetPrm(info->da,&prm);

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

1064:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
1065:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
1066:   MatSetOption(B,MAT_SYMMETRIC,PETSC_TRUE);
1067:   if (thi->verbose) {THIMatrixStatistics(thi,B,PETSC_VIEWER_STDOUT_WORLD);}
1068:   return(0);
1069: }

1073: static PetscErrorCode THIJacobianLocal_3D(DMDALocalInfo *info,Node ***x,Mat B,THI thi,THIAssemblyMode amode)
1074: {
1075:   PetscInt       xs,ys,xm,ym,zm,i,j,k,q,l,ll;
1076:   PetscReal      hx,hy;
1077:   PrmNode        **prm;

1081:   xs = info->zs;
1082:   ys = info->ys;
1083:   xm = info->zm;
1084:   ym = info->ym;
1085:   zm = info->xm;
1086:   hx = thi->Lx / info->mz;
1087:   hy = thi->Ly / info->my;

1089:   MatZeroEntries(B);
1090:   THIDAGetPrm(info->da,&prm);

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

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

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

1263:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
1264:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
1265:   MatSetOption(B,MAT_SYMMETRIC,PETSC_TRUE);
1266:   if (thi->verbose) {THIMatrixStatistics(thi,B,PETSC_VIEWER_STDOUT_WORLD);}
1267:   return(0);
1268: }

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

1277:   THIJacobianLocal_3D(info,x,B,thi,THIASSEMBLY_FULL);
1278:   return(0);
1279: }

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

1288:   THIJacobianLocal_3D(info,x,B,thi,THIASSEMBLY_TRIDIAGONAL);
1289:   return(0);
1290: }

1294: static PetscErrorCode DMRefineHierarchy_THI(DM dac0,PetscInt nlevels,DM hierarchy[])
1295: {
1296:   PetscErrorCode  ierr;
1297:   THI             thi;
1298:   PetscInt        dim,M,N,m,n,s,dof;
1299:   DM              dac,daf;
1300:   DMDAStencilType st;
1301:   DM_DA           *ddf,*ddc;

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

1315:   /* Creates a 3D DMDA with the same map-plane layout as the 2D one, with contiguous columns */
1316:   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);

1318:   daf->ops->creatematrix        = dac->ops->creatematrix;
1319:   daf->ops->createinterpolation = dac->ops->createinterpolation;
1320:   daf->ops->getcoloring         = dac->ops->getcoloring;
1321:   ddf                           = (DM_DA*)daf->data;
1322:   ddc                           = (DM_DA*)dac->data;
1323:   ddf->interptype               = ddc->interptype;

1325:   DMDASetFieldName(daf,0,"x-velocity");
1326:   DMDASetFieldName(daf,1,"y-velocity");

1328:   hierarchy[nlevels-1] = daf;
1329:   return(0);
1330: }

1334: static PetscErrorCode DMCreateInterpolation_DA_THI(DM dac,DM daf,Mat *A,Vec *scale)
1335: {
1337:   PetscInt       dim;

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

1352:     DMDAGetInfo(daf,0, &mz,&my,&mx, 0,0,0, 0,0,0,0,0,0);
1353:     DMDAGetCorners(daf,&zs,&ys,&xs,&zm,&ym,&xm);
1354:     if (zs != 0) SETERRQ(PETSC_COMM_SELF,1,"unexpected");
1355:     MatCreate(PetscObjectComm((PetscObject)daf),&B);
1356:     MatSetSizes(B,xm*ym*zm,xm*ym,mx*my*mz,mx*my);

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

1382: static PetscErrorCode DMCreateMatrix_THI_Tridiagonal(DM da,Mat *J)
1383: {
1384:   PetscErrorCode         ierr;
1385:   Mat                    A;
1386:   PetscInt               xm,ym,zm,dim,dof = 2,starts[3],dims[3];
1387:   ISLocalToGlobalMapping ltog;

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

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

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

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

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

1474:       PetscViewerASCIIPrintf(viewer,"      <PointData>\n");
1475:       PetscViewerASCIIPrintf(viewer,"        <DataArray type=\"Float32\" Name=\"velocity\" NumberOfComponents=\"3\" format=\"ascii\">\n");
1476:       for (i=0; i<nn; i+=dof) {
1477:         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);
1478:       }
1479:       PetscViewerASCIIPrintf(viewer,"        </DataArray>\n");

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

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

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

1512:   PetscInitialize(&argc,&argv,0,help);
1513:   comm = PETSC_COMM_WORLD;

1515:   THICreate(comm,&thi);
1516:   {
1517:     PetscInt M = 3,N = 3,P = 2;
1518:     PetscOptionsBegin(comm,NULL,"Grid resolution options","");
1519:     {
1520:       PetscOptionsInt("-M","Number of elements in x-direction on coarse level","",M,&M,NULL);
1521:       N    = M;
1522:       PetscOptionsInt("-N","Number of elements in y-direction on coarse level (if different from M)","",N,&N,NULL);
1523:       if (thi->coarse2d) {
1524:         PetscOptionsInt("-zlevels","Number of elements in z-direction on fine level","",thi->zlevels,&thi->zlevels,NULL);
1525:       } else {
1526:         PetscOptionsInt("-P","Number of elements in z-direction on coarse level","",P,&P,NULL);
1527:       }
1528:     }
1529:     PetscOptionsEnd();
1530:     if (thi->coarse2d) {
1531:       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);

1533:       da->ops->refinehierarchy     = DMRefineHierarchy_THI;
1534:       da->ops->createinterpolation = DMCreateInterpolation_DA_THI;

1536:       PetscObjectCompose((PetscObject)da,"THI",(PetscObject)thi);
1537:     } else {
1538:       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);
1539:     }
1540:     DMDASetFieldName(da,0,"x-velocity");
1541:     DMDASetFieldName(da,1,"y-velocity");
1542:   }
1543:   THISetUpDM(thi,da);
1544:   if (thi->tridiagonal) da->ops->creatematrix = DMCreateMatrix_THI_Tridiagonal;

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

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

1562:   DMSetApplicationContext(da,thi);

1564:   SNESCreate(comm,&snes);
1565:   SNESSetDM(snes,da);
1566:   DMDestroy(&da);
1567:   SNESSetComputeInitialGuess(snes,THIInitial,NULL);
1568:   SNESSetFromOptions(snes);

1570:   SNESSolve(snes,NULL,NULL);

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

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

1587:   DMDestroy(&da);
1588:   SNESDestroy(&snes);
1589:   THIDestroy(&thi);
1590:   PetscFinalize();
1591:   return 0;
1592: }