1*c4762a1bSJed Brown static const char help[] = "Toy hydrostatic ice flow with multigrid in 3D.\n\ 2*c4762a1bSJed Brown \n\ 3*c4762a1bSJed Brown Solves the hydrostatic (aka Blatter/Pattyn/First Order) equations for ice sheet flow\n\ 4*c4762a1bSJed Brown using multigrid. The ice uses a power-law rheology with \"Glen\" exponent 3 (corresponds\n\ 5*c4762a1bSJed Brown to p=4/3 in a p-Laplacian). The focus is on ISMIP-HOM experiments which assume periodic\n\ 6*c4762a1bSJed Brown boundary conditions in the x- and y-directions.\n\ 7*c4762a1bSJed Brown \n\ 8*c4762a1bSJed Brown Equations are rescaled so that the domain size and solution are O(1), details of this scaling\n\ 9*c4762a1bSJed Brown can be controlled by the options -units_meter, -units_second, and -units_kilogram.\n\ 10*c4762a1bSJed Brown \n\ 11*c4762a1bSJed Brown A VTK StructuredGrid output file can be written using the option -o filename.vts\n\ 12*c4762a1bSJed Brown \n\n"; 13*c4762a1bSJed Brown 14*c4762a1bSJed Brown /* 15*c4762a1bSJed Brown The equations for horizontal velocity (u,v) are 16*c4762a1bSJed Brown 17*c4762a1bSJed Brown - [eta (4 u_x + 2 v_y)]_x - [eta (u_y + v_x)]_y - [eta u_z]_z + rho g s_x = 0 18*c4762a1bSJed Brown - [eta (4 v_y + 2 u_x)]_y - [eta (u_y + v_x)]_x - [eta v_z]_z + rho g s_y = 0 19*c4762a1bSJed Brown 20*c4762a1bSJed Brown where 21*c4762a1bSJed Brown 22*c4762a1bSJed Brown eta = B/2 (epsilon + gamma)^((p-2)/2) 23*c4762a1bSJed Brown 24*c4762a1bSJed Brown is the nonlinear effective viscosity with regularization epsilon and hardness parameter B, 25*c4762a1bSJed Brown written in terms of the second invariant 26*c4762a1bSJed Brown 27*c4762a1bSJed Brown 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 28*c4762a1bSJed Brown 29*c4762a1bSJed Brown The surface boundary conditions are the natural conditions. The basal boundary conditions 30*c4762a1bSJed Brown are either no-slip, or Navier (linear) slip with spatially variant friction coefficient beta^2. 31*c4762a1bSJed Brown 32*c4762a1bSJed Brown In the code, the equations for (u,v) are multiplied through by 1/(rho g) so that residuals are O(1). 33*c4762a1bSJed Brown 34*c4762a1bSJed Brown The discretization is Q1 finite elements, managed by a DMDA. The grid is never distorted in the 35*c4762a1bSJed Brown map (x,y) plane, but the bed and surface may be bumpy. This is handled as usual in FEM, through 36*c4762a1bSJed Brown the Jacobian of the coordinate transformation from a reference element to the physical element. 37*c4762a1bSJed Brown 38*c4762a1bSJed Brown Since ice-flow is tightly coupled in the z-direction (within columns), the DMDA is managed 39*c4762a1bSJed Brown specially so that columns are never distributed, and are always contiguous in memory. 40*c4762a1bSJed Brown This amounts to reversing the meaning of X,Y,Z compared to the DMDA's internal interpretation, 41*c4762a1bSJed Brown and then indexing as vec[i][j][k]. The exotic coarse spaces require 2D DMDAs which are made to 42*c4762a1bSJed Brown use compatible domain decomposition relative to the 3D DMDAs. 43*c4762a1bSJed Brown 44*c4762a1bSJed Brown */ 45*c4762a1bSJed Brown 46*c4762a1bSJed Brown #include <petscts.h> 47*c4762a1bSJed Brown #include <petscdm.h> 48*c4762a1bSJed Brown #include <petscdmda.h> 49*c4762a1bSJed Brown #include <petscdmcomposite.h> 50*c4762a1bSJed Brown #include <ctype.h> /* toupper() */ 51*c4762a1bSJed Brown #include <petsc/private/petscimpl.h> 52*c4762a1bSJed Brown 53*c4762a1bSJed Brown #if defined __SSE2__ 54*c4762a1bSJed Brown # include <emmintrin.h> 55*c4762a1bSJed Brown #endif 56*c4762a1bSJed Brown 57*c4762a1bSJed Brown /* The SSE2 kernels are only for PetscScalar=double on architectures that support it */ 58*c4762a1bSJed Brown #define USE_SSE2_KERNELS (!defined NO_SSE2 \ 59*c4762a1bSJed Brown && !defined PETSC_USE_COMPLEX \ 60*c4762a1bSJed Brown && !defined PETSC_USE_REAL_SINGLE \ 61*c4762a1bSJed Brown && defined __SSE2__) 62*c4762a1bSJed Brown 63*c4762a1bSJed Brown #if !defined __STDC_VERSION__ || __STDC_VERSION__ < 199901L 64*c4762a1bSJed Brown # if defined __cplusplus /* C++ restrict is nonstandard and compilers have inconsistent rules about where it can be used */ 65*c4762a1bSJed Brown # define restrict 66*c4762a1bSJed Brown # else 67*c4762a1bSJed Brown # define restrict PETSC_RESTRICT 68*c4762a1bSJed Brown # endif 69*c4762a1bSJed Brown #endif 70*c4762a1bSJed Brown 71*c4762a1bSJed Brown static PetscClassId THI_CLASSID; 72*c4762a1bSJed Brown 73*c4762a1bSJed Brown typedef enum {QUAD_GAUSS,QUAD_LOBATTO} QuadratureType; 74*c4762a1bSJed Brown static const char *QuadratureTypes[] = {"gauss","lobatto","QuadratureType","QUAD_",0}; 75*c4762a1bSJed Brown static const PetscReal HexQWeights[8] = {1,1,1,1,1,1,1,1}; 76*c4762a1bSJed Brown static const PetscReal HexQNodes[] = {-0.57735026918962573, 0.57735026918962573}; 77*c4762a1bSJed Brown #define G 0.57735026918962573 78*c4762a1bSJed Brown #define H (0.5*(1.+G)) 79*c4762a1bSJed Brown #define L (0.5*(1.-G)) 80*c4762a1bSJed Brown #define M (-0.5) 81*c4762a1bSJed Brown #define P (0.5) 82*c4762a1bSJed Brown /* Special quadrature: Lobatto in horizontal, Gauss in vertical */ 83*c4762a1bSJed Brown static const PetscReal HexQInterp_Lobatto[8][8] = {{H,0,0,0,L,0,0,0}, 84*c4762a1bSJed Brown {0,H,0,0,0,L,0,0}, 85*c4762a1bSJed Brown {0,0,H,0,0,0,L,0}, 86*c4762a1bSJed Brown {0,0,0,H,0,0,0,L}, 87*c4762a1bSJed Brown {L,0,0,0,H,0,0,0}, 88*c4762a1bSJed Brown {0,L,0,0,0,H,0,0}, 89*c4762a1bSJed Brown {0,0,L,0,0,0,H,0}, 90*c4762a1bSJed Brown {0,0,0,L,0,0,0,H}}; 91*c4762a1bSJed Brown static const PetscReal HexQDeriv_Lobatto[8][8][3] = { 92*c4762a1bSJed Brown {{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} }, 93*c4762a1bSJed Brown {{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} }, 94*c4762a1bSJed Brown {{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} }, 95*c4762a1bSJed Brown {{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}}, 96*c4762a1bSJed Brown {{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} }, 97*c4762a1bSJed Brown {{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} }, 98*c4762a1bSJed Brown {{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} }, 99*c4762a1bSJed Brown {{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}}}; 100*c4762a1bSJed Brown /* Stanndard Gauss */ 101*c4762a1bSJed Brown 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}, 102*c4762a1bSJed Brown {L*H*H,H*H*H,H*L*H,L*L*H, L*H*L,H*H*L,H*L*L,L*L*L}, 103*c4762a1bSJed Brown {L*L*H,H*L*H,H*H*H,L*H*H, L*L*L,H*L*L,H*H*L,L*H*L}, 104*c4762a1bSJed Brown {H*L*H,L*L*H,L*H*H,H*H*H, H*L*L,L*L*L,L*H*L,H*H*L}, 105*c4762a1bSJed Brown {H*H*L,L*H*L,L*L*L,H*L*L, H*H*H,L*H*H,L*L*H,H*L*H}, 106*c4762a1bSJed Brown {L*H*L,H*H*L,H*L*L,L*L*L, L*H*H,H*H*H,H*L*H,L*L*H}, 107*c4762a1bSJed Brown {L*L*L,H*L*L,H*H*L,L*H*L, L*L*H,H*L*H,H*H*H,L*H*H}, 108*c4762a1bSJed Brown {H*L*L,L*L*L,L*H*L,H*H*L, H*L*H,L*L*H,L*H*H,H*H*H}}; 109*c4762a1bSJed Brown static const PetscReal HexQDeriv_Gauss[8][8][3] = { 110*c4762a1bSJed Brown {{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}}, 111*c4762a1bSJed Brown {{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}}, 112*c4762a1bSJed Brown {{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}}, 113*c4762a1bSJed Brown {{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}}, 114*c4762a1bSJed Brown {{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}}, 115*c4762a1bSJed Brown {{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}}, 116*c4762a1bSJed Brown {{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}}, 117*c4762a1bSJed Brown {{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}}}; 118*c4762a1bSJed Brown static const PetscReal (*HexQInterp)[8],(*HexQDeriv)[8][3]; 119*c4762a1bSJed Brown /* Standard 2x2 Gauss quadrature for the bottom layer. */ 120*c4762a1bSJed Brown static const PetscReal QuadQInterp[4][4] = {{H*H,L*H,L*L,H*L}, 121*c4762a1bSJed Brown {L*H,H*H,H*L,L*L}, 122*c4762a1bSJed Brown {L*L,H*L,H*H,L*H}, 123*c4762a1bSJed Brown {H*L,L*L,L*H,H*H}}; 124*c4762a1bSJed Brown static const PetscReal QuadQDeriv[4][4][2] = { 125*c4762a1bSJed Brown {{M*H,M*H},{P*H,M*L},{P*L,P*L},{M*L,P*H}}, 126*c4762a1bSJed Brown {{M*H,M*L},{P*H,M*H},{P*L,P*H},{M*L,P*L}}, 127*c4762a1bSJed Brown {{M*L,M*L},{P*L,M*H},{P*H,P*H},{M*H,P*L}}, 128*c4762a1bSJed Brown {{M*L,M*H},{P*L,M*L},{P*H,P*L},{M*H,P*H}}}; 129*c4762a1bSJed Brown #undef G 130*c4762a1bSJed Brown #undef H 131*c4762a1bSJed Brown #undef L 132*c4762a1bSJed Brown #undef M 133*c4762a1bSJed Brown #undef P 134*c4762a1bSJed Brown 135*c4762a1bSJed Brown #define HexExtract(x,i,j,k,n) do { \ 136*c4762a1bSJed Brown (n)[0] = (x)[i][j][k]; \ 137*c4762a1bSJed Brown (n)[1] = (x)[i+1][j][k]; \ 138*c4762a1bSJed Brown (n)[2] = (x)[i+1][j+1][k]; \ 139*c4762a1bSJed Brown (n)[3] = (x)[i][j+1][k]; \ 140*c4762a1bSJed Brown (n)[4] = (x)[i][j][k+1]; \ 141*c4762a1bSJed Brown (n)[5] = (x)[i+1][j][k+1]; \ 142*c4762a1bSJed Brown (n)[6] = (x)[i+1][j+1][k+1]; \ 143*c4762a1bSJed Brown (n)[7] = (x)[i][j+1][k+1]; \ 144*c4762a1bSJed Brown } while (0) 145*c4762a1bSJed Brown 146*c4762a1bSJed Brown #define HexExtractRef(x,i,j,k,n) do { \ 147*c4762a1bSJed Brown (n)[0] = &(x)[i][j][k]; \ 148*c4762a1bSJed Brown (n)[1] = &(x)[i+1][j][k]; \ 149*c4762a1bSJed Brown (n)[2] = &(x)[i+1][j+1][k]; \ 150*c4762a1bSJed Brown (n)[3] = &(x)[i][j+1][k]; \ 151*c4762a1bSJed Brown (n)[4] = &(x)[i][j][k+1]; \ 152*c4762a1bSJed Brown (n)[5] = &(x)[i+1][j][k+1]; \ 153*c4762a1bSJed Brown (n)[6] = &(x)[i+1][j+1][k+1]; \ 154*c4762a1bSJed Brown (n)[7] = &(x)[i][j+1][k+1]; \ 155*c4762a1bSJed Brown } while (0) 156*c4762a1bSJed Brown 157*c4762a1bSJed Brown #define QuadExtract(x,i,j,n) do { \ 158*c4762a1bSJed Brown (n)[0] = (x)[i][j]; \ 159*c4762a1bSJed Brown (n)[1] = (x)[i+1][j]; \ 160*c4762a1bSJed Brown (n)[2] = (x)[i+1][j+1]; \ 161*c4762a1bSJed Brown (n)[3] = (x)[i][j+1]; \ 162*c4762a1bSJed Brown } while (0) 163*c4762a1bSJed Brown 164*c4762a1bSJed Brown static PetscScalar Sqr(PetscScalar a) {return a*a;} 165*c4762a1bSJed Brown 166*c4762a1bSJed Brown static void HexGrad(const PetscReal dphi[][3],const PetscReal zn[],PetscReal dz[]) 167*c4762a1bSJed Brown { 168*c4762a1bSJed Brown PetscInt i; 169*c4762a1bSJed Brown dz[0] = dz[1] = dz[2] = 0; 170*c4762a1bSJed Brown for (i=0; i<8; i++) { 171*c4762a1bSJed Brown dz[0] += dphi[i][0] * zn[i]; 172*c4762a1bSJed Brown dz[1] += dphi[i][1] * zn[i]; 173*c4762a1bSJed Brown dz[2] += dphi[i][2] * zn[i]; 174*c4762a1bSJed Brown } 175*c4762a1bSJed Brown } 176*c4762a1bSJed Brown 177*c4762a1bSJed Brown static void HexComputeGeometry(PetscInt q,PetscReal hx,PetscReal hy,const PetscReal dz[restrict],PetscReal phi[restrict],PetscReal dphi[restrict][3],PetscReal *restrict jw) 178*c4762a1bSJed Brown { 179*c4762a1bSJed Brown const PetscReal 180*c4762a1bSJed Brown jac[3][3] = {{hx/2,0,0}, {0,hy/2,0}, {dz[0],dz[1],dz[2]}} 181*c4762a1bSJed Brown ,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]}} 182*c4762a1bSJed Brown ,jdet = jac[0][0]*jac[1][1]*jac[2][2]; 183*c4762a1bSJed Brown PetscInt i; 184*c4762a1bSJed Brown 185*c4762a1bSJed Brown for (i=0; i<8; i++) { 186*c4762a1bSJed Brown const PetscReal *dphir = HexQDeriv[q][i]; 187*c4762a1bSJed Brown phi[i] = HexQInterp[q][i]; 188*c4762a1bSJed Brown dphi[i][0] = dphir[0]*ijac[0][0] + dphir[1]*ijac[1][0] + dphir[2]*ijac[2][0]; 189*c4762a1bSJed Brown dphi[i][1] = dphir[0]*ijac[0][1] + dphir[1]*ijac[1][1] + dphir[2]*ijac[2][1]; 190*c4762a1bSJed Brown dphi[i][2] = dphir[0]*ijac[0][2] + dphir[1]*ijac[1][2] + dphir[2]*ijac[2][2]; 191*c4762a1bSJed Brown } 192*c4762a1bSJed Brown *jw = 1.0 * jdet; 193*c4762a1bSJed Brown } 194*c4762a1bSJed Brown 195*c4762a1bSJed Brown typedef struct _p_THI *THI; 196*c4762a1bSJed Brown typedef struct _n_Units *Units; 197*c4762a1bSJed Brown 198*c4762a1bSJed Brown typedef struct { 199*c4762a1bSJed Brown PetscScalar u,v; 200*c4762a1bSJed Brown } Node; 201*c4762a1bSJed Brown 202*c4762a1bSJed Brown typedef struct { 203*c4762a1bSJed Brown PetscScalar b; /* bed */ 204*c4762a1bSJed Brown PetscScalar h; /* thickness */ 205*c4762a1bSJed Brown PetscScalar beta2; /* friction */ 206*c4762a1bSJed Brown } PrmNode; 207*c4762a1bSJed Brown 208*c4762a1bSJed Brown #define FieldSize(ntype) ((PetscInt)(sizeof(ntype)/sizeof(PetscScalar))) 209*c4762a1bSJed Brown #define FieldOffset(ntype,member) ((PetscInt)(offsetof(ntype,member)/sizeof(PetscScalar))) 210*c4762a1bSJed Brown #define FieldIndex(ntype,i,member) ((PetscInt)((i)*FieldSize(ntype) + FieldOffset(ntype,member))) 211*c4762a1bSJed Brown #define NODE_SIZE FieldSize(Node) 212*c4762a1bSJed Brown #define PRMNODE_SIZE FieldSize(PrmNode) 213*c4762a1bSJed Brown 214*c4762a1bSJed Brown typedef struct { 215*c4762a1bSJed Brown PetscReal min,max,cmin,cmax; 216*c4762a1bSJed Brown } PRange; 217*c4762a1bSJed Brown 218*c4762a1bSJed Brown struct _p_THI { 219*c4762a1bSJed Brown PETSCHEADER(int); 220*c4762a1bSJed Brown void (*initialize)(THI,PetscReal x,PetscReal y,PrmNode *p); 221*c4762a1bSJed Brown PetscInt nlevels; 222*c4762a1bSJed Brown PetscInt zlevels; 223*c4762a1bSJed Brown PetscReal Lx,Ly,Lz; /* Model domain */ 224*c4762a1bSJed Brown PetscReal alpha; /* Bed angle */ 225*c4762a1bSJed Brown Units units; 226*c4762a1bSJed Brown PetscReal dirichlet_scale; 227*c4762a1bSJed Brown PetscReal ssa_friction_scale; 228*c4762a1bSJed Brown PetscReal inertia; 229*c4762a1bSJed Brown PRange eta; 230*c4762a1bSJed Brown PRange beta2; 231*c4762a1bSJed Brown struct { 232*c4762a1bSJed Brown PetscReal Bd2,eps,exponent,glen_n; 233*c4762a1bSJed Brown } viscosity; 234*c4762a1bSJed Brown struct { 235*c4762a1bSJed Brown PetscReal irefgam,eps2,exponent; 236*c4762a1bSJed Brown } friction; 237*c4762a1bSJed Brown struct { 238*c4762a1bSJed Brown PetscReal rate,exponent,refvel; 239*c4762a1bSJed Brown } erosion; 240*c4762a1bSJed Brown PetscReal rhog; 241*c4762a1bSJed Brown PetscBool no_slip; 242*c4762a1bSJed Brown PetscBool verbose; 243*c4762a1bSJed Brown char *mattype; 244*c4762a1bSJed Brown char *monitor_basename; 245*c4762a1bSJed Brown PetscInt monitor_interval; 246*c4762a1bSJed Brown }; 247*c4762a1bSJed Brown 248*c4762a1bSJed Brown struct _n_Units { 249*c4762a1bSJed Brown /* fundamental */ 250*c4762a1bSJed Brown PetscReal meter; 251*c4762a1bSJed Brown PetscReal kilogram; 252*c4762a1bSJed Brown PetscReal second; 253*c4762a1bSJed Brown /* derived */ 254*c4762a1bSJed Brown PetscReal Pascal; 255*c4762a1bSJed Brown PetscReal year; 256*c4762a1bSJed Brown }; 257*c4762a1bSJed Brown 258*c4762a1bSJed Brown static void PrmHexGetZ(const PrmNode pn[],PetscInt k,PetscInt zm,PetscReal zn[]) 259*c4762a1bSJed Brown { 260*c4762a1bSJed Brown const PetscScalar zm1 = zm-1, 261*c4762a1bSJed Brown znl[8] = {pn[0].b + pn[0].h*(PetscScalar)k/zm1, 262*c4762a1bSJed Brown pn[1].b + pn[1].h*(PetscScalar)k/zm1, 263*c4762a1bSJed Brown pn[2].b + pn[2].h*(PetscScalar)k/zm1, 264*c4762a1bSJed Brown pn[3].b + pn[3].h*(PetscScalar)k/zm1, 265*c4762a1bSJed Brown pn[0].b + pn[0].h*(PetscScalar)(k+1)/zm1, 266*c4762a1bSJed Brown pn[1].b + pn[1].h*(PetscScalar)(k+1)/zm1, 267*c4762a1bSJed Brown pn[2].b + pn[2].h*(PetscScalar)(k+1)/zm1, 268*c4762a1bSJed Brown pn[3].b + pn[3].h*(PetscScalar)(k+1)/zm1}; 269*c4762a1bSJed Brown PetscInt i; 270*c4762a1bSJed Brown for (i=0; i<8; i++) zn[i] = PetscRealPart(znl[i]); 271*c4762a1bSJed Brown } 272*c4762a1bSJed Brown 273*c4762a1bSJed Brown /* Compute a gradient of all the 2D fields at four quadrature points. Output for [quadrature_point][direction].field_name */ 274*c4762a1bSJed Brown static PetscErrorCode QuadComputeGrad4(const PetscReal dphi[][4][2],PetscReal hx,PetscReal hy,const PrmNode pn[4],PrmNode dp[4][2]) 275*c4762a1bSJed Brown { 276*c4762a1bSJed Brown PetscErrorCode ierr; 277*c4762a1bSJed Brown PetscInt q,i,f; 278*c4762a1bSJed Brown const PetscScalar (*restrict pg)[PRMNODE_SIZE] = (const PetscScalar(*)[PRMNODE_SIZE])pn; /* Get generic array pointers to the node */ 279*c4762a1bSJed Brown PetscScalar (*restrict dpg)[2][PRMNODE_SIZE] = (PetscScalar(*)[2][PRMNODE_SIZE])dp; 280*c4762a1bSJed Brown 281*c4762a1bSJed Brown PetscFunctionBeginUser; 282*c4762a1bSJed Brown ierr = PetscArrayzero(dpg,4);CHKERRQ(ierr); 283*c4762a1bSJed Brown for (q=0; q<4; q++) { 284*c4762a1bSJed Brown for (i=0; i<4; i++) { 285*c4762a1bSJed Brown for (f=0; f<PRMNODE_SIZE; f++) { 286*c4762a1bSJed Brown dpg[q][0][f] += dphi[q][i][0]/hx * pg[i][f]; 287*c4762a1bSJed Brown dpg[q][1][f] += dphi[q][i][1]/hy * pg[i][f]; 288*c4762a1bSJed Brown } 289*c4762a1bSJed Brown } 290*c4762a1bSJed Brown } 291*c4762a1bSJed Brown PetscFunctionReturn(0); 292*c4762a1bSJed Brown } 293*c4762a1bSJed Brown 294*c4762a1bSJed Brown static inline PetscReal StaggeredMidpoint2D(PetscScalar a,PetscScalar b,PetscScalar c,PetscScalar d) 295*c4762a1bSJed Brown {return 0.5*PetscRealPart(0.75*a + 0.75*b + 0.25*c + 0.25*d);} 296*c4762a1bSJed Brown static inline PetscReal UpwindFlux1D(PetscReal u,PetscReal hL,PetscReal hR) 297*c4762a1bSJed Brown {return (u > 0) ? hL*u : hR*u;} 298*c4762a1bSJed Brown 299*c4762a1bSJed Brown #define UpwindFluxXW(x3,x2,h,i,j,k,dj) UpwindFlux1D(StaggeredMidpoint2D(x3[i][j][k].u,x3[i-1][j][k].u, x3[i-1][j+dj][k].u,x3[i][k+dj][k].u), \ 300*c4762a1bSJed Brown PetscRealPart(0.75*x2[i-1][j ].h+0.25*x2[i-1][j+dj].h), PetscRealPart(0.75*x2[i ][j ].h+0.25*x2[i ][j+dj].h)) 301*c4762a1bSJed Brown #define UpwindFluxXE(x3,x2,h,i,j,k,dj) UpwindFlux1D(StaggeredMidpoint2D(x3[i][j][k].u,x3[i+1][j][k].u, x3[i+1][j+dj][k].u,x3[i][k+dj][k].u), \ 302*c4762a1bSJed Brown PetscRealPart(0.75*x2[i ][j ].h+0.25*x2[i ][j+dj].h), PetscRealPart(0.75*x2[i+1][j ].h+0.25*x2[i+1][j+dj].h)) 303*c4762a1bSJed Brown #define UpwindFluxYS(x3,x2,h,i,j,k,di) UpwindFlux1D(StaggeredMidpoint2D(x3[i][j][k].v,x3[i][j-1][k].v, x3[i+di][j-1][k].v,x3[i+di][j][k].v), \ 304*c4762a1bSJed Brown PetscRealPart(0.75*x2[i ][j-1].h+0.25*x2[i+di][j-1].h), PetscRealPart(0.75*x2[i ][j ].h+0.25*x2[i+di][j ].h)) 305*c4762a1bSJed Brown #define UpwindFluxYN(x3,x2,h,i,j,k,di) UpwindFlux1D(StaggeredMidpoint2D(x3[i][j][k].v,x3[i][j+1][k].v, x3[i+di][j+1][k].v,x3[i+di][j][k].v), \ 306*c4762a1bSJed Brown PetscRealPart(0.75*x2[i ][j ].h+0.25*x2[i+di][j ].h), PetscRealPart(0.75*x2[i ][j+1].h+0.25*x2[i+di][j+1].h)) 307*c4762a1bSJed Brown 308*c4762a1bSJed Brown static void PrmNodeGetFaceMeasure(const PrmNode **p,PetscInt i,PetscInt j,PetscScalar h[]) 309*c4762a1bSJed Brown { 310*c4762a1bSJed Brown /* West */ 311*c4762a1bSJed Brown h[0] = StaggeredMidpoint2D(p[i][j].h,p[i-1][j].h,p[i-1][j-1].h,p[i][j-1].h); 312*c4762a1bSJed Brown h[1] = StaggeredMidpoint2D(p[i][j].h,p[i-1][j].h,p[i-1][j+1].h,p[i][j+1].h); 313*c4762a1bSJed Brown /* East */ 314*c4762a1bSJed Brown h[2] = StaggeredMidpoint2D(p[i][j].h,p[i+1][j].h,p[i+1][j+1].h,p[i][j+1].h); 315*c4762a1bSJed Brown h[3] = StaggeredMidpoint2D(p[i][j].h,p[i+1][j].h,p[i+1][j-1].h,p[i][j-1].h); 316*c4762a1bSJed Brown /* South */ 317*c4762a1bSJed Brown h[4] = StaggeredMidpoint2D(p[i][j].h,p[i][j-1].h,p[i+1][j-1].h,p[i+1][j].h); 318*c4762a1bSJed Brown h[5] = StaggeredMidpoint2D(p[i][j].h,p[i][j-1].h,p[i-1][j-1].h,p[i-1][j].h); 319*c4762a1bSJed Brown /* North */ 320*c4762a1bSJed Brown h[6] = StaggeredMidpoint2D(p[i][j].h,p[i][j+1].h,p[i-1][j+1].h,p[i-1][j].h); 321*c4762a1bSJed Brown h[7] = StaggeredMidpoint2D(p[i][j].h,p[i][j+1].h,p[i+1][j+1].h,p[i+1][j].h); 322*c4762a1bSJed Brown } 323*c4762a1bSJed Brown 324*c4762a1bSJed Brown /* Tests A and C are from the ISMIP-HOM paper (Pattyn et al. 2008) */ 325*c4762a1bSJed Brown static void THIInitialize_HOM_A(THI thi,PetscReal x,PetscReal y,PrmNode *p) 326*c4762a1bSJed Brown { 327*c4762a1bSJed Brown Units units = thi->units; 328*c4762a1bSJed Brown PetscReal s = -x*PetscSinReal(thi->alpha); 329*c4762a1bSJed Brown p->b = s - 1000*units->meter + 500*units->meter * PetscSinReal(x*2*PETSC_PI/thi->Lx) * PetscSinReal(y*2*PETSC_PI/thi->Ly); 330*c4762a1bSJed Brown p->h = s - p->b; 331*c4762a1bSJed Brown p->beta2 = -1e-10; /* This value is not used, but it should not be huge because that would change the finite difference step size */ 332*c4762a1bSJed Brown } 333*c4762a1bSJed Brown 334*c4762a1bSJed Brown static void THIInitialize_HOM_C(THI thi,PetscReal x,PetscReal y,PrmNode *p) 335*c4762a1bSJed Brown { 336*c4762a1bSJed Brown Units units = thi->units; 337*c4762a1bSJed Brown PetscReal s = -x*PetscSinReal(thi->alpha); 338*c4762a1bSJed Brown p->b = s - 1000*units->meter; 339*c4762a1bSJed Brown p->h = s - p->b; 340*c4762a1bSJed Brown /* tau_b = beta2 v is a stress (Pa). 341*c4762a1bSJed Brown * This is a big number in our units (it needs to balance the driving force from the surface), so we scale it by 1/rhog, just like the residual. */ 342*c4762a1bSJed Brown p->beta2 = 1000 * (1 + PetscSinReal(x*2*PETSC_PI/thi->Lx)*PetscSinReal(y*2*PETSC_PI/thi->Ly)) * units->Pascal * units->year / units->meter / thi->rhog; 343*c4762a1bSJed Brown } 344*c4762a1bSJed Brown 345*c4762a1bSJed Brown /* These are just toys */ 346*c4762a1bSJed Brown 347*c4762a1bSJed Brown /* From Fred Herman */ 348*c4762a1bSJed Brown static void THIInitialize_HOM_F(THI thi,PetscReal x,PetscReal y,PrmNode *p) 349*c4762a1bSJed Brown { 350*c4762a1bSJed Brown Units units = thi->units; 351*c4762a1bSJed Brown PetscReal s = -x*PetscSinReal(thi->alpha); 352*c4762a1bSJed Brown p->b = s - 1000*units->meter + 100*units->meter * PetscSinReal(x*2*PETSC_PI/thi->Lx);/* * sin(y*2*PETSC_PI/thi->Ly); */ 353*c4762a1bSJed Brown p->h = s - p->b; 354*c4762a1bSJed Brown p->h = (1-(atan((x-thi->Lx/2)/1.)+PETSC_PI/2.)/PETSC_PI)*500*units->meter+1*units->meter; 355*c4762a1bSJed Brown s = PetscRealPart(p->b + p->h); 356*c4762a1bSJed Brown p->beta2 = -1e-10; 357*c4762a1bSJed Brown /* p->beta2 = 1000 * units->Pascal * units->year / units->meter; */ 358*c4762a1bSJed Brown } 359*c4762a1bSJed Brown 360*c4762a1bSJed Brown /* Same bed as test A, free slip everywhere except for a discontinuous jump to a circular sticky region in the middle. */ 361*c4762a1bSJed Brown static void THIInitialize_HOM_X(THI thi,PetscReal xx,PetscReal yy,PrmNode *p) 362*c4762a1bSJed Brown { 363*c4762a1bSJed Brown Units units = thi->units; 364*c4762a1bSJed Brown PetscReal x = xx*2*PETSC_PI/thi->Lx - PETSC_PI,y = yy*2*PETSC_PI/thi->Ly - PETSC_PI; /* [-pi,pi] */ 365*c4762a1bSJed Brown PetscReal r = PetscSqrtReal(x*x + y*y),s = -x*PetscSinReal(thi->alpha); 366*c4762a1bSJed Brown p->b = s - 1000*units->meter + 500*units->meter * PetscSinReal(x + PETSC_PI) * PetscSinReal(y + PETSC_PI); 367*c4762a1bSJed Brown p->h = s - p->b; 368*c4762a1bSJed Brown p->beta2 = 1000 * (r < 1 ? 2 : 0) * units->Pascal * units->year / units->meter / thi->rhog; 369*c4762a1bSJed Brown } 370*c4762a1bSJed Brown 371*c4762a1bSJed Brown /* Like Z, but with 200 meter cliffs */ 372*c4762a1bSJed Brown static void THIInitialize_HOM_Y(THI thi,PetscReal xx,PetscReal yy,PrmNode *p) 373*c4762a1bSJed Brown { 374*c4762a1bSJed Brown Units units = thi->units; 375*c4762a1bSJed Brown PetscReal x = xx*2*PETSC_PI/thi->Lx - PETSC_PI,y = yy*2*PETSC_PI/thi->Ly - PETSC_PI; /* [-pi,pi] */ 376*c4762a1bSJed Brown PetscReal r = PetscSqrtReal(x*x + y*y),s = -x*PetscSinReal(thi->alpha); 377*c4762a1bSJed Brown p->b = s - 1000*units->meter + 500*units->meter * PetscSinReal(x + PETSC_PI) * PetscSinReal(y + PETSC_PI); 378*c4762a1bSJed Brown if (PetscRealPart(p->b) > -700*units->meter) p->b += 200*units->meter; 379*c4762a1bSJed Brown p->h = s - p->b; 380*c4762a1bSJed Brown 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 / thi->rhog; 381*c4762a1bSJed Brown } 382*c4762a1bSJed Brown 383*c4762a1bSJed Brown /* Same bed as A, smoothly varying slipperiness, similar to MATLAB's "sombrero" (uncorrelated with bathymetry) */ 384*c4762a1bSJed Brown static void THIInitialize_HOM_Z(THI thi,PetscReal xx,PetscReal yy,PrmNode *p) 385*c4762a1bSJed Brown { 386*c4762a1bSJed Brown Units units = thi->units; 387*c4762a1bSJed Brown PetscReal x = xx*2*PETSC_PI/thi->Lx - PETSC_PI,y = yy*2*PETSC_PI/thi->Ly - PETSC_PI; /* [-pi,pi] */ 388*c4762a1bSJed Brown PetscReal r = PetscSqrtReal(x*x + y*y),s = -x*PetscSinReal(thi->alpha); 389*c4762a1bSJed Brown p->b = s - 1000*units->meter + 500*units->meter * PetscSinReal(x + PETSC_PI) * PetscSinReal(y + PETSC_PI); 390*c4762a1bSJed Brown p->h = s - p->b; 391*c4762a1bSJed Brown 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 / thi->rhog; 392*c4762a1bSJed Brown } 393*c4762a1bSJed Brown 394*c4762a1bSJed Brown static void THIFriction(THI thi,PetscReal rbeta2,PetscReal gam,PetscReal *beta2,PetscReal *dbeta2) 395*c4762a1bSJed Brown { 396*c4762a1bSJed Brown if (thi->friction.irefgam == 0) { 397*c4762a1bSJed Brown Units units = thi->units; 398*c4762a1bSJed Brown thi->friction.irefgam = 1./(0.5*PetscSqr(100 * units->meter / units->year)); 399*c4762a1bSJed Brown thi->friction.eps2 = 0.5*PetscSqr(1.e-4 / thi->friction.irefgam); 400*c4762a1bSJed Brown } 401*c4762a1bSJed Brown if (thi->friction.exponent == 0) { 402*c4762a1bSJed Brown *beta2 = rbeta2; 403*c4762a1bSJed Brown *dbeta2 = 0; 404*c4762a1bSJed Brown } else { 405*c4762a1bSJed Brown *beta2 = rbeta2 * PetscPowReal(thi->friction.eps2 + gam*thi->friction.irefgam,thi->friction.exponent); 406*c4762a1bSJed Brown *dbeta2 = thi->friction.exponent * *beta2 / (thi->friction.eps2 + gam*thi->friction.irefgam) * thi->friction.irefgam; 407*c4762a1bSJed Brown } 408*c4762a1bSJed Brown } 409*c4762a1bSJed Brown 410*c4762a1bSJed Brown static void THIViscosity(THI thi,PetscReal gam,PetscReal *eta,PetscReal *deta) 411*c4762a1bSJed Brown { 412*c4762a1bSJed Brown PetscReal Bd2,eps,exponent; 413*c4762a1bSJed Brown if (thi->viscosity.Bd2 == 0) { 414*c4762a1bSJed Brown Units units = thi->units; 415*c4762a1bSJed Brown const PetscReal 416*c4762a1bSJed Brown n = thi->viscosity.glen_n, /* Glen exponent */ 417*c4762a1bSJed Brown p = 1. + 1./n, /* for Stokes */ 418*c4762a1bSJed Brown A = 1.e-16 * PetscPowReal(units->Pascal,-n) / units->year, /* softness parameter (Pa^{-n}/s) */ 419*c4762a1bSJed Brown B = PetscPowReal(A,-1./n); /* hardness parameter */ 420*c4762a1bSJed Brown thi->viscosity.Bd2 = B/2; 421*c4762a1bSJed Brown thi->viscosity.exponent = (p-2)/2; 422*c4762a1bSJed Brown thi->viscosity.eps = 0.5*PetscSqr(1e-5 / units->year); 423*c4762a1bSJed Brown } 424*c4762a1bSJed Brown Bd2 = thi->viscosity.Bd2; 425*c4762a1bSJed Brown exponent = thi->viscosity.exponent; 426*c4762a1bSJed Brown eps = thi->viscosity.eps; 427*c4762a1bSJed Brown *eta = Bd2 * PetscPowReal(eps + gam,exponent); 428*c4762a1bSJed Brown *deta = exponent * (*eta) / (eps + gam); 429*c4762a1bSJed Brown } 430*c4762a1bSJed Brown 431*c4762a1bSJed Brown static void THIErosion(THI thi,const Node *vel,PetscScalar *erate,Node *derate) 432*c4762a1bSJed Brown { 433*c4762a1bSJed Brown const PetscScalar magref2 = 1.e-10 + (PetscSqr(vel->u) + PetscSqr(vel->v)) / PetscSqr(thi->erosion.refvel), 434*c4762a1bSJed Brown rate = -thi->erosion.rate*PetscPowScalar(magref2, 0.5*thi->erosion.exponent); 435*c4762a1bSJed Brown if (erate) *erate = rate; 436*c4762a1bSJed Brown if (derate) { 437*c4762a1bSJed Brown if (thi->erosion.exponent == 1) { 438*c4762a1bSJed Brown derate->u = 0; 439*c4762a1bSJed Brown derate->v = 0; 440*c4762a1bSJed Brown } else { 441*c4762a1bSJed Brown derate->u = 0.5*thi->erosion.exponent * rate / magref2 * 2. * vel->u / PetscSqr(thi->erosion.refvel); 442*c4762a1bSJed Brown derate->v = 0.5*thi->erosion.exponent * rate / magref2 * 2. * vel->v / PetscSqr(thi->erosion.refvel); 443*c4762a1bSJed Brown } 444*c4762a1bSJed Brown } 445*c4762a1bSJed Brown } 446*c4762a1bSJed Brown 447*c4762a1bSJed Brown static void RangeUpdate(PetscReal *min,PetscReal *max,PetscReal x) 448*c4762a1bSJed Brown { 449*c4762a1bSJed Brown if (x < *min) *min = x; 450*c4762a1bSJed Brown if (x > *max) *max = x; 451*c4762a1bSJed Brown } 452*c4762a1bSJed Brown 453*c4762a1bSJed Brown static void PRangeClear(PRange *p) 454*c4762a1bSJed Brown { 455*c4762a1bSJed Brown p->cmin = p->min = 1e100; 456*c4762a1bSJed Brown p->cmax = p->max = -1e100; 457*c4762a1bSJed Brown } 458*c4762a1bSJed Brown 459*c4762a1bSJed Brown static PetscErrorCode PRangeMinMax(PRange *p,PetscReal min,PetscReal max) 460*c4762a1bSJed Brown { 461*c4762a1bSJed Brown 462*c4762a1bSJed Brown PetscFunctionBeginUser; 463*c4762a1bSJed Brown p->cmin = min; 464*c4762a1bSJed Brown p->cmax = max; 465*c4762a1bSJed Brown if (min < p->min) p->min = min; 466*c4762a1bSJed Brown if (max > p->max) p->max = max; 467*c4762a1bSJed Brown PetscFunctionReturn(0); 468*c4762a1bSJed Brown } 469*c4762a1bSJed Brown 470*c4762a1bSJed Brown static PetscErrorCode THIDestroy(THI *thi) 471*c4762a1bSJed Brown { 472*c4762a1bSJed Brown PetscErrorCode ierr; 473*c4762a1bSJed Brown 474*c4762a1bSJed Brown PetscFunctionBeginUser; 475*c4762a1bSJed Brown if (--((PetscObject)(*thi))->refct > 0) PetscFunctionReturn(0); 476*c4762a1bSJed Brown ierr = PetscFree((*thi)->units);CHKERRQ(ierr); 477*c4762a1bSJed Brown ierr = PetscFree((*thi)->mattype);CHKERRQ(ierr); 478*c4762a1bSJed Brown ierr = PetscFree((*thi)->monitor_basename);CHKERRQ(ierr); 479*c4762a1bSJed Brown ierr = PetscHeaderDestroy(thi);CHKERRQ(ierr); 480*c4762a1bSJed Brown PetscFunctionReturn(0); 481*c4762a1bSJed Brown } 482*c4762a1bSJed Brown 483*c4762a1bSJed Brown static PetscErrorCode THICreate(MPI_Comm comm,THI *inthi) 484*c4762a1bSJed Brown { 485*c4762a1bSJed Brown static PetscBool registered = PETSC_FALSE; 486*c4762a1bSJed Brown THI thi; 487*c4762a1bSJed Brown Units units; 488*c4762a1bSJed Brown char monitor_basename[PETSC_MAX_PATH_LEN] = "thi-"; 489*c4762a1bSJed Brown PetscErrorCode ierr; 490*c4762a1bSJed Brown 491*c4762a1bSJed Brown PetscFunctionBeginUser; 492*c4762a1bSJed Brown *inthi = 0; 493*c4762a1bSJed Brown if (!registered) { 494*c4762a1bSJed Brown ierr = PetscClassIdRegister("Toy Hydrostatic Ice",&THI_CLASSID);CHKERRQ(ierr); 495*c4762a1bSJed Brown registered = PETSC_TRUE; 496*c4762a1bSJed Brown } 497*c4762a1bSJed Brown ierr = PetscHeaderCreate(thi,THI_CLASSID,"THI","Toy Hydrostatic Ice","THI",comm,THIDestroy,0);CHKERRQ(ierr); 498*c4762a1bSJed Brown 499*c4762a1bSJed Brown ierr = PetscNew(&thi->units);CHKERRQ(ierr); 500*c4762a1bSJed Brown 501*c4762a1bSJed Brown units = thi->units; 502*c4762a1bSJed Brown units->meter = 1e-2; 503*c4762a1bSJed Brown units->second = 1e-7; 504*c4762a1bSJed Brown units->kilogram = 1e-12; 505*c4762a1bSJed Brown 506*c4762a1bSJed Brown ierr = PetscOptionsBegin(comm,NULL,"Scaled units options","");CHKERRQ(ierr); 507*c4762a1bSJed Brown { 508*c4762a1bSJed Brown ierr = PetscOptionsReal("-units_meter","1 meter in scaled length units","",units->meter,&units->meter,NULL);CHKERRQ(ierr); 509*c4762a1bSJed Brown ierr = PetscOptionsReal("-units_second","1 second in scaled time units","",units->second,&units->second,NULL);CHKERRQ(ierr); 510*c4762a1bSJed Brown ierr = PetscOptionsReal("-units_kilogram","1 kilogram in scaled mass units","",units->kilogram,&units->kilogram,NULL);CHKERRQ(ierr); 511*c4762a1bSJed Brown } 512*c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 513*c4762a1bSJed Brown units->Pascal = units->kilogram / (units->meter * PetscSqr(units->second)); 514*c4762a1bSJed Brown units->year = 31556926. * units->second, /* seconds per year */ 515*c4762a1bSJed Brown 516*c4762a1bSJed Brown thi->Lx = 10.e3; 517*c4762a1bSJed Brown thi->Ly = 10.e3; 518*c4762a1bSJed Brown thi->Lz = 1000; 519*c4762a1bSJed Brown thi->nlevels = 1; 520*c4762a1bSJed Brown thi->dirichlet_scale = 1; 521*c4762a1bSJed Brown thi->verbose = PETSC_FALSE; 522*c4762a1bSJed Brown 523*c4762a1bSJed Brown thi->viscosity.glen_n = 3.; 524*c4762a1bSJed Brown thi->erosion.rate = 1e-3; /* m/a */ 525*c4762a1bSJed Brown thi->erosion.exponent = 1.; 526*c4762a1bSJed Brown thi->erosion.refvel = 1.; /* m/a */ 527*c4762a1bSJed Brown 528*c4762a1bSJed Brown ierr = PetscOptionsBegin(comm,NULL,"Toy Hydrostatic Ice options","");CHKERRQ(ierr); 529*c4762a1bSJed Brown { 530*c4762a1bSJed Brown QuadratureType quad = QUAD_GAUSS; 531*c4762a1bSJed Brown char homexp[] = "A"; 532*c4762a1bSJed Brown char mtype[256] = MATSBAIJ; 533*c4762a1bSJed Brown PetscReal L,m = 1.0; 534*c4762a1bSJed Brown PetscBool flg; 535*c4762a1bSJed Brown L = thi->Lx; 536*c4762a1bSJed Brown ierr = PetscOptionsReal("-thi_L","Domain size (m)","",L,&L,&flg);CHKERRQ(ierr); 537*c4762a1bSJed Brown if (flg) thi->Lx = thi->Ly = L; 538*c4762a1bSJed Brown ierr = PetscOptionsReal("-thi_Lx","X Domain size (m)","",thi->Lx,&thi->Lx,NULL);CHKERRQ(ierr); 539*c4762a1bSJed Brown ierr = PetscOptionsReal("-thi_Ly","Y Domain size (m)","",thi->Ly,&thi->Ly,NULL);CHKERRQ(ierr); 540*c4762a1bSJed Brown ierr = PetscOptionsReal("-thi_Lz","Z Domain size (m)","",thi->Lz,&thi->Lz,NULL);CHKERRQ(ierr); 541*c4762a1bSJed Brown ierr = PetscOptionsString("-thi_hom","ISMIP-HOM experiment (A or C)","",homexp,homexp,sizeof(homexp),NULL);CHKERRQ(ierr); 542*c4762a1bSJed Brown switch (homexp[0] = toupper(homexp[0])) { 543*c4762a1bSJed Brown case 'A': 544*c4762a1bSJed Brown thi->initialize = THIInitialize_HOM_A; 545*c4762a1bSJed Brown thi->no_slip = PETSC_TRUE; 546*c4762a1bSJed Brown thi->alpha = 0.5; 547*c4762a1bSJed Brown break; 548*c4762a1bSJed Brown case 'C': 549*c4762a1bSJed Brown thi->initialize = THIInitialize_HOM_C; 550*c4762a1bSJed Brown thi->no_slip = PETSC_FALSE; 551*c4762a1bSJed Brown thi->alpha = 0.1; 552*c4762a1bSJed Brown break; 553*c4762a1bSJed Brown case 'F': 554*c4762a1bSJed Brown thi->initialize = THIInitialize_HOM_F; 555*c4762a1bSJed Brown thi->no_slip = PETSC_FALSE; 556*c4762a1bSJed Brown thi->alpha = 0.5; 557*c4762a1bSJed Brown break; 558*c4762a1bSJed Brown case 'X': 559*c4762a1bSJed Brown thi->initialize = THIInitialize_HOM_X; 560*c4762a1bSJed Brown thi->no_slip = PETSC_FALSE; 561*c4762a1bSJed Brown thi->alpha = 0.3; 562*c4762a1bSJed Brown break; 563*c4762a1bSJed Brown case 'Y': 564*c4762a1bSJed Brown thi->initialize = THIInitialize_HOM_Y; 565*c4762a1bSJed Brown thi->no_slip = PETSC_FALSE; 566*c4762a1bSJed Brown thi->alpha = 0.5; 567*c4762a1bSJed Brown break; 568*c4762a1bSJed Brown case 'Z': 569*c4762a1bSJed Brown thi->initialize = THIInitialize_HOM_Z; 570*c4762a1bSJed Brown thi->no_slip = PETSC_FALSE; 571*c4762a1bSJed Brown thi->alpha = 0.5; 572*c4762a1bSJed Brown break; 573*c4762a1bSJed Brown default: 574*c4762a1bSJed Brown SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"HOM experiment '%c' not implemented",homexp[0]); 575*c4762a1bSJed Brown } 576*c4762a1bSJed Brown ierr = PetscOptionsEnum("-thi_quadrature","Quadrature to use for 3D elements","",QuadratureTypes,(PetscEnum)quad,(PetscEnum*)&quad,NULL);CHKERRQ(ierr); 577*c4762a1bSJed Brown switch (quad) { 578*c4762a1bSJed Brown case QUAD_GAUSS: 579*c4762a1bSJed Brown HexQInterp = HexQInterp_Gauss; 580*c4762a1bSJed Brown HexQDeriv = HexQDeriv_Gauss; 581*c4762a1bSJed Brown break; 582*c4762a1bSJed Brown case QUAD_LOBATTO: 583*c4762a1bSJed Brown HexQInterp = HexQInterp_Lobatto; 584*c4762a1bSJed Brown HexQDeriv = HexQDeriv_Lobatto; 585*c4762a1bSJed Brown break; 586*c4762a1bSJed Brown } 587*c4762a1bSJed Brown ierr = PetscOptionsReal("-thi_alpha","Bed angle (degrees)","",thi->alpha,&thi->alpha,NULL);CHKERRQ(ierr); 588*c4762a1bSJed Brown ierr = PetscOptionsReal("-thi_viscosity_glen_n","Exponent in Glen flow law, 1=linear, infty=ideal plastic",NULL,thi->viscosity.glen_n,&thi->viscosity.glen_n,NULL);CHKERRQ(ierr); 589*c4762a1bSJed Brown ierr = PetscOptionsReal("-thi_friction_m","Friction exponent, 0=Coulomb, 1=Navier","",m,&m,NULL);CHKERRQ(ierr); 590*c4762a1bSJed Brown thi->friction.exponent = (m-1)/2; 591*c4762a1bSJed Brown ierr = PetscOptionsReal("-thi_erosion_rate","Rate of erosion relative to sliding velocity at reference velocity (m/a)",NULL,thi->erosion.rate,&thi->erosion.rate,NULL);CHKERRQ(ierr); 592*c4762a1bSJed Brown ierr = PetscOptionsReal("-thi_erosion_exponent","Power of sliding velocity appearing in erosion relation",NULL,thi->erosion.exponent,&thi->erosion.exponent,NULL);CHKERRQ(ierr); 593*c4762a1bSJed Brown ierr = PetscOptionsReal("-thi_erosion_refvel","Reference sliding velocity for erosion (m/a)",NULL,thi->erosion.refvel,&thi->erosion.refvel,NULL);CHKERRQ(ierr); 594*c4762a1bSJed Brown thi->erosion.rate *= units->meter / units->year; 595*c4762a1bSJed Brown thi->erosion.refvel *= units->meter / units->year; 596*c4762a1bSJed Brown ierr = PetscOptionsReal("-thi_dirichlet_scale","Scale Dirichlet boundary conditions by this factor","",thi->dirichlet_scale,&thi->dirichlet_scale,NULL);CHKERRQ(ierr); 597*c4762a1bSJed Brown ierr = 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);CHKERRQ(ierr); 598*c4762a1bSJed Brown ierr = PetscOptionsReal("-thi_inertia","Coefficient of accelaration term in velocity system, physical is almost zero",NULL,thi->inertia,&thi->inertia,NULL);CHKERRQ(ierr); 599*c4762a1bSJed Brown ierr = PetscOptionsInt("-thi_nlevels","Number of levels of refinement","",thi->nlevels,&thi->nlevels,NULL);CHKERRQ(ierr); 600*c4762a1bSJed Brown ierr = PetscOptionsFList("-thi_mat_type","Matrix type","MatSetType",MatList,mtype,(char*)mtype,sizeof(mtype),NULL);CHKERRQ(ierr); 601*c4762a1bSJed Brown ierr = PetscStrallocpy(mtype,&thi->mattype);CHKERRQ(ierr); 602*c4762a1bSJed Brown ierr = PetscOptionsBool("-thi_verbose","Enable verbose output (like matrix sizes and statistics)","",thi->verbose,&thi->verbose,NULL);CHKERRQ(ierr); 603*c4762a1bSJed Brown ierr = PetscOptionsString("-thi_monitor","Basename to write state files to",NULL,monitor_basename,monitor_basename,sizeof(monitor_basename),&flg);CHKERRQ(ierr); 604*c4762a1bSJed Brown if (flg) { 605*c4762a1bSJed Brown ierr = PetscStrallocpy(monitor_basename,&thi->monitor_basename);CHKERRQ(ierr); 606*c4762a1bSJed Brown thi->monitor_interval = 1; 607*c4762a1bSJed Brown ierr = PetscOptionsInt("-thi_monitor_interval","Frequency at which to write state files",NULL,thi->monitor_interval,&thi->monitor_interval,NULL);CHKERRQ(ierr); 608*c4762a1bSJed Brown } 609*c4762a1bSJed Brown } 610*c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 611*c4762a1bSJed Brown 612*c4762a1bSJed Brown /* dimensionalize */ 613*c4762a1bSJed Brown thi->Lx *= units->meter; 614*c4762a1bSJed Brown thi->Ly *= units->meter; 615*c4762a1bSJed Brown thi->Lz *= units->meter; 616*c4762a1bSJed Brown thi->alpha *= PETSC_PI / 180; 617*c4762a1bSJed Brown 618*c4762a1bSJed Brown PRangeClear(&thi->eta); 619*c4762a1bSJed Brown PRangeClear(&thi->beta2); 620*c4762a1bSJed Brown 621*c4762a1bSJed Brown { 622*c4762a1bSJed Brown PetscReal u = 1000*units->meter/(3e7*units->second), 623*c4762a1bSJed Brown gradu = u / (100*units->meter),eta,deta, 624*c4762a1bSJed Brown rho = 910 * units->kilogram/PetscPowRealInt(units->meter,3), 625*c4762a1bSJed Brown grav = 9.81 * units->meter/PetscSqr(units->second), 626*c4762a1bSJed Brown driving = rho * grav * PetscSinReal(thi->alpha) * 1000*units->meter; 627*c4762a1bSJed Brown THIViscosity(thi,0.5*gradu*gradu,&eta,&deta); 628*c4762a1bSJed Brown thi->rhog = rho * grav; 629*c4762a1bSJed Brown if (thi->verbose) { 630*c4762a1bSJed Brown ierr = PetscPrintf(PetscObjectComm((PetscObject)thi),"Units: meter %8.2g second %8.2g kg %8.2g Pa %8.2g\n",units->meter,units->second,units->kilogram,units->Pascal);CHKERRQ(ierr); 631*c4762a1bSJed Brown ierr = PetscPrintf(PetscObjectComm((PetscObject)thi),"Domain (%6.2g,%6.2g,%6.2g), pressure %8.2g, driving stress %8.2g\n",thi->Lx,thi->Ly,thi->Lz,rho*grav*1e3*units->meter,driving);CHKERRQ(ierr); 632*c4762a1bSJed Brown ierr = PetscPrintf(PetscObjectComm((PetscObject)thi),"Large velocity 1km/a %8.2g, velocity gradient %8.2g, eta %8.2g, stress %8.2g, ratio %8.2g\n",u,gradu,eta,2*eta*gradu,2*eta*gradu/driving);CHKERRQ(ierr); 633*c4762a1bSJed Brown THIViscosity(thi,0.5*PetscSqr(1e-3*gradu),&eta,&deta); 634*c4762a1bSJed Brown ierr = PetscPrintf(PetscObjectComm((PetscObject)thi),"Small velocity 1m/a %8.2g, velocity gradient %8.2g, eta %8.2g, stress %8.2g, ratio %8.2g\n",1e-3*u,1e-3*gradu,eta,2*eta*1e-3*gradu,2*eta*1e-3*gradu/driving);CHKERRQ(ierr); 635*c4762a1bSJed Brown } 636*c4762a1bSJed Brown } 637*c4762a1bSJed Brown 638*c4762a1bSJed Brown *inthi = thi; 639*c4762a1bSJed Brown PetscFunctionReturn(0); 640*c4762a1bSJed Brown } 641*c4762a1bSJed Brown 642*c4762a1bSJed Brown /* Our problem is periodic, but the domain has a mean slope of alpha so the bed does not line up between the upstream 643*c4762a1bSJed Brown * and downstream ends of the domain. This function fixes the ghost values so that the domain appears truly periodic in 644*c4762a1bSJed Brown * the horizontal. */ 645*c4762a1bSJed Brown static PetscErrorCode THIFixGhosts(THI thi,DM da3,DM da2,Vec X3,Vec X2) 646*c4762a1bSJed Brown { 647*c4762a1bSJed Brown PetscErrorCode ierr; 648*c4762a1bSJed Brown DMDALocalInfo info; 649*c4762a1bSJed Brown PrmNode **x2; 650*c4762a1bSJed Brown PetscInt i,j; 651*c4762a1bSJed Brown 652*c4762a1bSJed Brown PetscFunctionBeginUser; 653*c4762a1bSJed Brown ierr = DMDAGetLocalInfo(da3,&info);CHKERRQ(ierr); 654*c4762a1bSJed Brown /* ierr = VecView(X2,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */ 655*c4762a1bSJed Brown ierr = DMDAVecGetArray(da2,X2,&x2);CHKERRQ(ierr); 656*c4762a1bSJed Brown for (i=info.gzs; i<info.gzs+info.gzm; i++) { 657*c4762a1bSJed Brown if (i > -1 && i < info.mz) continue; 658*c4762a1bSJed Brown for (j=info.gys; j<info.gys+info.gym; j++) { 659*c4762a1bSJed Brown x2[i][j].b += PetscSinReal(thi->alpha) * thi->Lx * (i<0 ? 1.0 : -1.0); 660*c4762a1bSJed Brown } 661*c4762a1bSJed Brown } 662*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da2,X2,&x2);CHKERRQ(ierr); 663*c4762a1bSJed Brown /* ierr = VecView(X2,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */ 664*c4762a1bSJed Brown PetscFunctionReturn(0); 665*c4762a1bSJed Brown } 666*c4762a1bSJed Brown 667*c4762a1bSJed Brown static PetscErrorCode THIInitializePrm(THI thi,DM da2prm,PrmNode **p) 668*c4762a1bSJed Brown { 669*c4762a1bSJed Brown PetscInt i,j,xs,xm,ys,ym,mx,my; 670*c4762a1bSJed Brown PetscErrorCode ierr; 671*c4762a1bSJed Brown 672*c4762a1bSJed Brown PetscFunctionBeginUser; 673*c4762a1bSJed Brown ierr = DMDAGetGhostCorners(da2prm,&ys,&xs,0,&ym,&xm,0);CHKERRQ(ierr); 674*c4762a1bSJed Brown ierr = DMDAGetInfo(da2prm,0, &my,&mx,0, 0,0,0, 0,0,0,0,0,0);CHKERRQ(ierr); 675*c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 676*c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 677*c4762a1bSJed Brown PetscReal xx = thi->Lx*i/mx,yy = thi->Ly*j/my; 678*c4762a1bSJed Brown thi->initialize(thi,xx,yy,&p[i][j]); 679*c4762a1bSJed Brown } 680*c4762a1bSJed Brown } 681*c4762a1bSJed Brown PetscFunctionReturn(0); 682*c4762a1bSJed Brown } 683*c4762a1bSJed Brown 684*c4762a1bSJed Brown static PetscErrorCode THIInitial(THI thi,DM pack,Vec X) 685*c4762a1bSJed Brown { 686*c4762a1bSJed Brown DM da3,da2; 687*c4762a1bSJed Brown PetscInt i,j,k,xs,xm,ys,ym,zs,zm,mx,my; 688*c4762a1bSJed Brown PetscReal hx,hy; 689*c4762a1bSJed Brown PrmNode **prm; 690*c4762a1bSJed Brown Node ***x; 691*c4762a1bSJed Brown Vec X3g,X2g,X2; 692*c4762a1bSJed Brown PetscErrorCode ierr; 693*c4762a1bSJed Brown 694*c4762a1bSJed Brown PetscFunctionBeginUser; 695*c4762a1bSJed Brown ierr = DMCompositeGetEntries(pack,&da3,&da2);CHKERRQ(ierr); 696*c4762a1bSJed Brown ierr = DMCompositeGetAccess(pack,X,&X3g,&X2g);CHKERRQ(ierr); 697*c4762a1bSJed Brown ierr = DMGetLocalVector(da2,&X2);CHKERRQ(ierr); 698*c4762a1bSJed Brown 699*c4762a1bSJed Brown ierr = DMDAGetInfo(da3,0, 0,&my,&mx, 0,0,0, 0,0,0,0,0,0);CHKERRQ(ierr); 700*c4762a1bSJed Brown ierr = DMDAGetCorners(da3,&zs,&ys,&xs,&zm,&ym,&xm);CHKERRQ(ierr); 701*c4762a1bSJed Brown ierr = DMDAVecGetArray(da3,X3g,&x);CHKERRQ(ierr); 702*c4762a1bSJed Brown ierr = DMDAVecGetArray(da2,X2,&prm);CHKERRQ(ierr); 703*c4762a1bSJed Brown 704*c4762a1bSJed Brown ierr = THIInitializePrm(thi,da2,prm);CHKERRQ(ierr); 705*c4762a1bSJed Brown 706*c4762a1bSJed Brown hx = thi->Lx / mx; 707*c4762a1bSJed Brown hy = thi->Ly / my; 708*c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 709*c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 710*c4762a1bSJed Brown for (k=zs; k<zs+zm; k++) { 711*c4762a1bSJed Brown const PetscScalar zm1 = zm-1, 712*c4762a1bSJed Brown 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), 713*c4762a1bSJed Brown 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); 714*c4762a1bSJed Brown x[i][j][k].u = 0. * drivingx * prm[i][j].h*(PetscScalar)k/zm1; 715*c4762a1bSJed Brown x[i][j][k].v = 0. * drivingy * prm[i][j].h*(PetscScalar)k/zm1; 716*c4762a1bSJed Brown } 717*c4762a1bSJed Brown } 718*c4762a1bSJed Brown } 719*c4762a1bSJed Brown 720*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da3,X3g,&x);CHKERRQ(ierr); 721*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da2,X2,&prm);CHKERRQ(ierr); 722*c4762a1bSJed Brown 723*c4762a1bSJed Brown ierr = DMLocalToGlobalBegin(da2,X2,INSERT_VALUES,X2g);CHKERRQ(ierr); 724*c4762a1bSJed Brown ierr = DMLocalToGlobalEnd (da2,X2,INSERT_VALUES,X2g);CHKERRQ(ierr); 725*c4762a1bSJed Brown ierr = DMRestoreLocalVector(da2,&X2);CHKERRQ(ierr); 726*c4762a1bSJed Brown 727*c4762a1bSJed Brown ierr = DMCompositeRestoreAccess(pack,X,&X3g,&X2g);CHKERRQ(ierr); 728*c4762a1bSJed Brown PetscFunctionReturn(0); 729*c4762a1bSJed Brown } 730*c4762a1bSJed Brown 731*c4762a1bSJed Brown static void PointwiseNonlinearity(THI thi,const Node n[restrict 8],const PetscReal phi[restrict 3],PetscReal dphi[restrict 8][3],PetscScalar *restrict u,PetscScalar *restrict v,PetscScalar du[restrict 3],PetscScalar dv[restrict 3],PetscReal *eta,PetscReal *deta) 732*c4762a1bSJed Brown { 733*c4762a1bSJed Brown PetscInt l,ll; 734*c4762a1bSJed Brown PetscScalar gam; 735*c4762a1bSJed Brown 736*c4762a1bSJed Brown du[0] = du[1] = du[2] = 0; 737*c4762a1bSJed Brown dv[0] = dv[1] = dv[2] = 0; 738*c4762a1bSJed Brown *u = 0; 739*c4762a1bSJed Brown *v = 0; 740*c4762a1bSJed Brown for (l=0; l<8; l++) { 741*c4762a1bSJed Brown *u += phi[l] * n[l].u; 742*c4762a1bSJed Brown *v += phi[l] * n[l].v; 743*c4762a1bSJed Brown for (ll=0; ll<3; ll++) { 744*c4762a1bSJed Brown du[ll] += dphi[l][ll] * n[l].u; 745*c4762a1bSJed Brown dv[ll] += dphi[l][ll] * n[l].v; 746*c4762a1bSJed Brown } 747*c4762a1bSJed Brown } 748*c4762a1bSJed Brown 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]); 749*c4762a1bSJed Brown THIViscosity(thi,PetscRealPart(gam),eta,deta); 750*c4762a1bSJed Brown } 751*c4762a1bSJed Brown 752*c4762a1bSJed Brown static PetscErrorCode THIFunctionLocal_3D(DMDALocalInfo *info,const Node ***x,const PrmNode **prm,const Node ***xdot,Node ***f,THI thi) 753*c4762a1bSJed Brown { 754*c4762a1bSJed Brown PetscInt xs,ys,xm,ym,zm,i,j,k,q,l; 755*c4762a1bSJed Brown PetscReal hx,hy,etamin,etamax,beta2min,beta2max; 756*c4762a1bSJed Brown PetscErrorCode ierr; 757*c4762a1bSJed Brown 758*c4762a1bSJed Brown PetscFunctionBeginUser; 759*c4762a1bSJed Brown xs = info->zs; 760*c4762a1bSJed Brown ys = info->ys; 761*c4762a1bSJed Brown xm = info->zm; 762*c4762a1bSJed Brown ym = info->ym; 763*c4762a1bSJed Brown zm = info->xm; 764*c4762a1bSJed Brown hx = thi->Lx / info->mz; 765*c4762a1bSJed Brown hy = thi->Ly / info->my; 766*c4762a1bSJed Brown 767*c4762a1bSJed Brown etamin = 1e100; 768*c4762a1bSJed Brown etamax = 0; 769*c4762a1bSJed Brown beta2min = 1e100; 770*c4762a1bSJed Brown beta2max = 0; 771*c4762a1bSJed Brown 772*c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 773*c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 774*c4762a1bSJed Brown PrmNode pn[4],dpn[4][2]; 775*c4762a1bSJed Brown QuadExtract(prm,i,j,pn); 776*c4762a1bSJed Brown ierr = QuadComputeGrad4(QuadQDeriv,hx,hy,pn,dpn);CHKERRQ(ierr); 777*c4762a1bSJed Brown for (k=0; k<zm-1; k++) { 778*c4762a1bSJed Brown PetscInt ls = 0; 779*c4762a1bSJed Brown Node n[8],ndot[8],*fn[8]; 780*c4762a1bSJed Brown PetscReal zn[8],etabase = 0; 781*c4762a1bSJed Brown PrmHexGetZ(pn,k,zm,zn); 782*c4762a1bSJed Brown HexExtract(x,i,j,k,n); 783*c4762a1bSJed Brown HexExtract(xdot,i,j,k,ndot);CHKERRQ(ierr); 784*c4762a1bSJed Brown HexExtractRef(f,i,j,k,fn); 785*c4762a1bSJed Brown if (thi->no_slip && k == 0) { 786*c4762a1bSJed Brown for (l=0; l<4; l++) n[l].u = n[l].v = 0; 787*c4762a1bSJed Brown /* The first 4 basis functions lie on the bottom layer, so their contribution is exactly 0, hence we can skip them */ 788*c4762a1bSJed Brown ls = 4; 789*c4762a1bSJed Brown } 790*c4762a1bSJed Brown for (q=0; q<8; q++) { 791*c4762a1bSJed Brown PetscReal dz[3],phi[8],dphi[8][3],jw,eta,deta; 792*c4762a1bSJed Brown PetscScalar du[3],dv[3],u,v,udot=0,vdot=0; 793*c4762a1bSJed Brown for (l=ls; l<8; l++) { 794*c4762a1bSJed Brown udot += HexQInterp[q][l]*ndot[l].u; 795*c4762a1bSJed Brown vdot += HexQInterp[q][l]*ndot[l].v; 796*c4762a1bSJed Brown } 797*c4762a1bSJed Brown HexGrad(HexQDeriv[q],zn,dz); 798*c4762a1bSJed Brown HexComputeGeometry(q,hx,hy,dz,phi,dphi,&jw); 799*c4762a1bSJed Brown PointwiseNonlinearity(thi,n,phi,dphi,&u,&v,du,dv,&eta,&deta); 800*c4762a1bSJed Brown jw /= thi->rhog; /* scales residuals to be O(1) */ 801*c4762a1bSJed Brown if (q == 0) etabase = eta; 802*c4762a1bSJed Brown RangeUpdate(&etamin,&etamax,eta); 803*c4762a1bSJed Brown for (l=ls; l<8; l++) { /* test functions */ 804*c4762a1bSJed Brown const PetscScalar ds[2] = {dpn[q%4][0].h+dpn[q%4][0].b, dpn[q%4][1].h+dpn[q%4][1].b}; 805*c4762a1bSJed Brown const PetscReal pp = phi[l],*dp = dphi[l]; 806*c4762a1bSJed Brown 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]; 807*c4762a1bSJed Brown 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]; 808*c4762a1bSJed Brown fn[l]->u += pp*jw*udot*thi->inertia*pp; 809*c4762a1bSJed Brown fn[l]->v += pp*jw*vdot*thi->inertia*pp; 810*c4762a1bSJed Brown } 811*c4762a1bSJed Brown } 812*c4762a1bSJed Brown if (k == 0) { /* we are on a bottom face */ 813*c4762a1bSJed Brown if (thi->no_slip) { 814*c4762a1bSJed Brown /* Note: Non-Galerkin coarse grid operators are very sensitive to the scaling of Dirichlet boundary 815*c4762a1bSJed Brown * conditions. After shenanigans above, etabase contains the effective viscosity at the closest quadrature 816*c4762a1bSJed Brown * point to the bed. We want the diagonal entry in the Dirichlet condition to have similar magnitude to the 817*c4762a1bSJed Brown * diagonal entry corresponding to the adjacent node. The fundamental scaling of the viscous part is in 818*c4762a1bSJed Brown * diagu, diagv below. This scaling is easy to recognize by considering the finite difference operator after 819*c4762a1bSJed Brown * scaling by element size. The no-slip Dirichlet condition is scaled by this factor, and also in the 820*c4762a1bSJed Brown * assembled matrix (see the similar block in THIJacobianLocal). 821*c4762a1bSJed Brown * 822*c4762a1bSJed Brown * Note that the residual at this Dirichlet node is linear in the state at this node, but also depends 823*c4762a1bSJed Brown * (nonlinearly in general) on the neighboring interior nodes through the local viscosity. This will make 824*c4762a1bSJed Brown * a matrix-free Jacobian have extra entries in the corresponding row. We assemble only the diagonal part, 825*c4762a1bSJed Brown * so the solution will exactly satisfy the boundary condition after the first linear iteration. 826*c4762a1bSJed Brown */ 827*c4762a1bSJed Brown const PetscReal hz = PetscRealPart(pn[0].h)/(zm-1.); 828*c4762a1bSJed Brown 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); 829*c4762a1bSJed Brown fn[0]->u = thi->dirichlet_scale*diagu*x[i][j][k].u; 830*c4762a1bSJed Brown fn[0]->v = thi->dirichlet_scale*diagv*x[i][j][k].v; 831*c4762a1bSJed Brown } else { /* Integrate over bottom face to apply boundary condition */ 832*c4762a1bSJed Brown for (q=0; q<4; q++) { /* We remove the explicit scaling of the residual by 1/rhog because beta2 already has that scaling to be O(1) */ 833*c4762a1bSJed Brown const PetscReal jw = 0.25*hx*hy,*phi = QuadQInterp[q]; 834*c4762a1bSJed Brown PetscScalar u =0,v=0,rbeta2=0; 835*c4762a1bSJed Brown PetscReal beta2,dbeta2; 836*c4762a1bSJed Brown for (l=0; l<4; l++) { 837*c4762a1bSJed Brown u += phi[l]*n[l].u; 838*c4762a1bSJed Brown v += phi[l]*n[l].v; 839*c4762a1bSJed Brown rbeta2 += phi[l]*pn[l].beta2; 840*c4762a1bSJed Brown } 841*c4762a1bSJed Brown THIFriction(thi,PetscRealPart(rbeta2),PetscRealPart(u*u+v*v)/2,&beta2,&dbeta2); 842*c4762a1bSJed Brown RangeUpdate(&beta2min,&beta2max,beta2); 843*c4762a1bSJed Brown for (l=0; l<4; l++) { 844*c4762a1bSJed Brown const PetscReal pp = phi[l]; 845*c4762a1bSJed Brown fn[ls+l]->u += pp*jw*beta2*u; 846*c4762a1bSJed Brown fn[ls+l]->v += pp*jw*beta2*v; 847*c4762a1bSJed Brown } 848*c4762a1bSJed Brown } 849*c4762a1bSJed Brown } 850*c4762a1bSJed Brown } 851*c4762a1bSJed Brown } 852*c4762a1bSJed Brown } 853*c4762a1bSJed Brown } 854*c4762a1bSJed Brown 855*c4762a1bSJed Brown ierr = PRangeMinMax(&thi->eta,etamin,etamax);CHKERRQ(ierr); 856*c4762a1bSJed Brown ierr = PRangeMinMax(&thi->beta2,beta2min,beta2max);CHKERRQ(ierr); 857*c4762a1bSJed Brown PetscFunctionReturn(0); 858*c4762a1bSJed Brown } 859*c4762a1bSJed Brown 860*c4762a1bSJed Brown static PetscErrorCode THIFunctionLocal_2D(DMDALocalInfo *info,const Node ***x,const PrmNode **prm,const PrmNode **prmdot,PrmNode **f,THI thi) 861*c4762a1bSJed Brown { 862*c4762a1bSJed Brown PetscInt xs,ys,xm,ym,zm,i,j,k; 863*c4762a1bSJed Brown 864*c4762a1bSJed Brown PetscFunctionBeginUser; 865*c4762a1bSJed Brown xs = info->zs; 866*c4762a1bSJed Brown ys = info->ys; 867*c4762a1bSJed Brown xm = info->zm; 868*c4762a1bSJed Brown ym = info->ym; 869*c4762a1bSJed Brown zm = info->xm; 870*c4762a1bSJed Brown 871*c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 872*c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 873*c4762a1bSJed Brown PetscScalar div = 0,erate,h[8]; 874*c4762a1bSJed Brown PrmNodeGetFaceMeasure(prm,i,j,h); 875*c4762a1bSJed Brown for (k=0; k<zm; k++) { 876*c4762a1bSJed Brown PetscScalar weight = (k==0 || k == zm-1) ? 0.5/(zm-1) : 1.0/(zm-1); 877*c4762a1bSJed Brown if (0) { /* centered flux */ 878*c4762a1bSJed Brown div += (- weight*h[0] * StaggeredMidpoint2D(x[i][j][k].u,x[i-1][j][k].u, x[i-1][j-1][k].u,x[i][j-1][k].u) 879*c4762a1bSJed Brown - weight*h[1] * StaggeredMidpoint2D(x[i][j][k].u,x[i-1][j][k].u, x[i-1][j+1][k].u,x[i][j+1][k].u) 880*c4762a1bSJed Brown + weight*h[2] * StaggeredMidpoint2D(x[i][j][k].u,x[i+1][j][k].u, x[i+1][j+1][k].u,x[i][j+1][k].u) 881*c4762a1bSJed Brown + weight*h[3] * StaggeredMidpoint2D(x[i][j][k].u,x[i+1][j][k].u, x[i+1][j-1][k].u,x[i][j-1][k].u) 882*c4762a1bSJed Brown - weight*h[4] * StaggeredMidpoint2D(x[i][j][k].v,x[i][j-1][k].v, x[i+1][j-1][k].v,x[i+1][j][k].v) 883*c4762a1bSJed Brown - weight*h[5] * StaggeredMidpoint2D(x[i][j][k].v,x[i][j-1][k].v, x[i-1][j-1][k].v,x[i-1][j][k].v) 884*c4762a1bSJed Brown + weight*h[6] * StaggeredMidpoint2D(x[i][j][k].v,x[i][j+1][k].v, x[i-1][j+1][k].v,x[i-1][j][k].v) 885*c4762a1bSJed Brown + weight*h[7] * StaggeredMidpoint2D(x[i][j][k].v,x[i][j+1][k].v, x[i+1][j+1][k].v,x[i+1][j][k].v)); 886*c4762a1bSJed Brown } else { /* Upwind flux */ 887*c4762a1bSJed Brown div += weight*(-UpwindFluxXW(x,prm,h,i,j,k, 1) 888*c4762a1bSJed Brown -UpwindFluxXW(x,prm,h,i,j,k,-1) 889*c4762a1bSJed Brown +UpwindFluxXE(x,prm,h,i,j,k, 1) 890*c4762a1bSJed Brown +UpwindFluxXE(x,prm,h,i,j,k,-1) 891*c4762a1bSJed Brown -UpwindFluxYS(x,prm,h,i,j,k, 1) 892*c4762a1bSJed Brown -UpwindFluxYS(x,prm,h,i,j,k,-1) 893*c4762a1bSJed Brown +UpwindFluxYN(x,prm,h,i,j,k, 1) 894*c4762a1bSJed Brown +UpwindFluxYN(x,prm,h,i,j,k,-1)); 895*c4762a1bSJed Brown } 896*c4762a1bSJed Brown } 897*c4762a1bSJed Brown /* printf("div[%d][%d] %g\n",i,j,div); */ 898*c4762a1bSJed Brown THIErosion(thi,&x[i][j][0],&erate,NULL); 899*c4762a1bSJed Brown f[i][j].b = prmdot[i][j].b - erate; 900*c4762a1bSJed Brown f[i][j].h = prmdot[i][j].h + div; 901*c4762a1bSJed Brown f[i][j].beta2 = prmdot[i][j].beta2; 902*c4762a1bSJed Brown } 903*c4762a1bSJed Brown } 904*c4762a1bSJed Brown PetscFunctionReturn(0); 905*c4762a1bSJed Brown } 906*c4762a1bSJed Brown 907*c4762a1bSJed Brown static PetscErrorCode THIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx) 908*c4762a1bSJed Brown { 909*c4762a1bSJed Brown PetscErrorCode ierr; 910*c4762a1bSJed Brown THI thi = (THI)ctx; 911*c4762a1bSJed Brown DM pack,da3,da2; 912*c4762a1bSJed Brown Vec X3,X2,Xdot3,Xdot2,F3,F2,F3g,F2g; 913*c4762a1bSJed Brown const Node ***x3,***xdot3; 914*c4762a1bSJed Brown const PrmNode **x2,**xdot2; 915*c4762a1bSJed Brown Node ***f3; 916*c4762a1bSJed Brown PrmNode **f2; 917*c4762a1bSJed Brown DMDALocalInfo info3; 918*c4762a1bSJed Brown 919*c4762a1bSJed Brown PetscFunctionBeginUser; 920*c4762a1bSJed Brown ierr = TSGetDM(ts,&pack);CHKERRQ(ierr); 921*c4762a1bSJed Brown ierr = DMCompositeGetEntries(pack,&da3,&da2);CHKERRQ(ierr); 922*c4762a1bSJed Brown ierr = DMDAGetLocalInfo(da3,&info3);CHKERRQ(ierr); 923*c4762a1bSJed Brown ierr = DMCompositeGetLocalVectors(pack,&X3,&X2);CHKERRQ(ierr); 924*c4762a1bSJed Brown ierr = DMCompositeGetLocalVectors(pack,&Xdot3,&Xdot2);CHKERRQ(ierr); 925*c4762a1bSJed Brown ierr = DMCompositeScatter(pack,X,X3,X2);CHKERRQ(ierr); 926*c4762a1bSJed Brown ierr = THIFixGhosts(thi,da3,da2,X3,X2);CHKERRQ(ierr); 927*c4762a1bSJed Brown ierr = DMCompositeScatter(pack,Xdot,Xdot3,Xdot2);CHKERRQ(ierr); 928*c4762a1bSJed Brown 929*c4762a1bSJed Brown ierr = DMGetLocalVector(da3,&F3);CHKERRQ(ierr); 930*c4762a1bSJed Brown ierr = DMGetLocalVector(da2,&F2);CHKERRQ(ierr); 931*c4762a1bSJed Brown ierr = VecZeroEntries(F3);CHKERRQ(ierr); 932*c4762a1bSJed Brown 933*c4762a1bSJed Brown ierr = DMDAVecGetArray(da3,X3,&x3);CHKERRQ(ierr); 934*c4762a1bSJed Brown ierr = DMDAVecGetArray(da2,X2,&x2);CHKERRQ(ierr); 935*c4762a1bSJed Brown ierr = DMDAVecGetArray(da3,Xdot3,&xdot3);CHKERRQ(ierr); 936*c4762a1bSJed Brown ierr = DMDAVecGetArray(da2,Xdot2,&xdot2);CHKERRQ(ierr); 937*c4762a1bSJed Brown ierr = DMDAVecGetArray(da3,F3,&f3);CHKERRQ(ierr); 938*c4762a1bSJed Brown ierr = DMDAVecGetArray(da2,F2,&f2);CHKERRQ(ierr); 939*c4762a1bSJed Brown 940*c4762a1bSJed Brown ierr = THIFunctionLocal_3D(&info3,x3,x2,xdot3,f3,thi);CHKERRQ(ierr); 941*c4762a1bSJed Brown ierr = THIFunctionLocal_2D(&info3,x3,x2,xdot2,f2,thi);CHKERRQ(ierr); 942*c4762a1bSJed Brown 943*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da3,X3,&x3);CHKERRQ(ierr); 944*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da2,X2,&x2);CHKERRQ(ierr); 945*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da3,Xdot3,&xdot3);CHKERRQ(ierr); 946*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da2,Xdot2,&xdot2);CHKERRQ(ierr); 947*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da3,F3,&f3);CHKERRQ(ierr); 948*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da2,F2,&f2);CHKERRQ(ierr); 949*c4762a1bSJed Brown 950*c4762a1bSJed Brown ierr = DMCompositeRestoreLocalVectors(pack,&X3,&X2);CHKERRQ(ierr); 951*c4762a1bSJed Brown ierr = DMCompositeRestoreLocalVectors(pack,&Xdot3,&Xdot2);CHKERRQ(ierr); 952*c4762a1bSJed Brown 953*c4762a1bSJed Brown ierr = VecZeroEntries(F);CHKERRQ(ierr); 954*c4762a1bSJed Brown ierr = DMCompositeGetAccess(pack,F,&F3g,&F2g);CHKERRQ(ierr); 955*c4762a1bSJed Brown ierr = DMLocalToGlobalBegin(da3,F3,ADD_VALUES,F3g);CHKERRQ(ierr); 956*c4762a1bSJed Brown ierr = DMLocalToGlobalEnd (da3,F3,ADD_VALUES,F3g);CHKERRQ(ierr); 957*c4762a1bSJed Brown ierr = DMLocalToGlobalBegin(da2,F2,INSERT_VALUES,F2g);CHKERRQ(ierr); 958*c4762a1bSJed Brown ierr = DMLocalToGlobalEnd (da2,F2,INSERT_VALUES,F2g);CHKERRQ(ierr); 959*c4762a1bSJed Brown 960*c4762a1bSJed Brown if (thi->verbose) { 961*c4762a1bSJed Brown PetscViewer viewer; 962*c4762a1bSJed Brown ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)thi),&viewer);CHKERRQ(ierr); 963*c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer,"3D_Velocity residual (bs=2):\n");CHKERRQ(ierr); 964*c4762a1bSJed Brown ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 965*c4762a1bSJed Brown ierr = VecView(F3,viewer);CHKERRQ(ierr); 966*c4762a1bSJed Brown ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 967*c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer,"2D_Fields residual (bs=3):\n");CHKERRQ(ierr); 968*c4762a1bSJed Brown ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 969*c4762a1bSJed Brown ierr = VecView(F2,viewer);CHKERRQ(ierr); 970*c4762a1bSJed Brown ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 971*c4762a1bSJed Brown } 972*c4762a1bSJed Brown 973*c4762a1bSJed Brown ierr = DMCompositeRestoreAccess(pack,F,&F3g,&F2g);CHKERRQ(ierr); 974*c4762a1bSJed Brown 975*c4762a1bSJed Brown ierr = DMRestoreLocalVector(da3,&F3);CHKERRQ(ierr); 976*c4762a1bSJed Brown ierr = DMRestoreLocalVector(da2,&F2);CHKERRQ(ierr); 977*c4762a1bSJed Brown PetscFunctionReturn(0); 978*c4762a1bSJed Brown } 979*c4762a1bSJed Brown 980*c4762a1bSJed Brown static PetscErrorCode THIMatrixStatistics(THI thi,Mat B,PetscViewer viewer) 981*c4762a1bSJed Brown { 982*c4762a1bSJed Brown PetscErrorCode ierr; 983*c4762a1bSJed Brown PetscReal nrm; 984*c4762a1bSJed Brown PetscInt m; 985*c4762a1bSJed Brown PetscMPIInt rank; 986*c4762a1bSJed Brown 987*c4762a1bSJed Brown PetscFunctionBeginUser; 988*c4762a1bSJed Brown ierr = MatNorm(B,NORM_FROBENIUS,&nrm);CHKERRQ(ierr); 989*c4762a1bSJed Brown ierr = MatGetSize(B,&m,0);CHKERRQ(ierr); 990*c4762a1bSJed Brown ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)B),&rank);CHKERRQ(ierr); 991*c4762a1bSJed Brown if (!rank) { 992*c4762a1bSJed Brown PetscScalar val0,val2; 993*c4762a1bSJed Brown ierr = MatGetValue(B,0,0,&val0);CHKERRQ(ierr); 994*c4762a1bSJed Brown ierr = MatGetValue(B,2,2,&val2);CHKERRQ(ierr); 995*c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer,"Matrix dim %8d norm %8.2e, (0,0) %8.2e (2,2) %8.2e, eta [%8.2e,%8.2e] beta2 [%8.2e,%8.2e]\n",m,nrm,PetscRealPart(val0),PetscRealPart(val2),thi->eta.cmin,thi->eta.cmax,thi->beta2.cmin,thi->beta2.cmax);CHKERRQ(ierr); 996*c4762a1bSJed Brown } 997*c4762a1bSJed Brown PetscFunctionReturn(0); 998*c4762a1bSJed Brown } 999*c4762a1bSJed Brown 1000*c4762a1bSJed Brown static PetscErrorCode THISurfaceStatistics(DM pack,Vec X,PetscReal *min,PetscReal *max,PetscReal *mean) 1001*c4762a1bSJed Brown { 1002*c4762a1bSJed Brown PetscErrorCode ierr; 1003*c4762a1bSJed Brown DM da3,da2; 1004*c4762a1bSJed Brown Vec X3,X2; 1005*c4762a1bSJed Brown Node ***x; 1006*c4762a1bSJed Brown PetscInt i,j,xs,ys,zs,xm,ym,zm,mx,my,mz; 1007*c4762a1bSJed Brown PetscReal umin = 1e100,umax=-1e100; 1008*c4762a1bSJed Brown PetscScalar usum =0.0,gusum; 1009*c4762a1bSJed Brown 1010*c4762a1bSJed Brown PetscFunctionBeginUser; 1011*c4762a1bSJed Brown ierr = DMCompositeGetEntries(pack,&da3,&da2);CHKERRQ(ierr); 1012*c4762a1bSJed Brown ierr = DMCompositeGetAccess(pack,X,&X3,&X2);CHKERRQ(ierr); 1013*c4762a1bSJed Brown *min = *max = *mean = 0; 1014*c4762a1bSJed Brown ierr = DMDAGetInfo(da3,0, &mz,&my,&mx, 0,0,0, 0,0,0,0,0,0);CHKERRQ(ierr); 1015*c4762a1bSJed Brown ierr = DMDAGetCorners(da3,&zs,&ys,&xs,&zm,&ym,&xm);CHKERRQ(ierr); 1016*c4762a1bSJed Brown if (zs != 0 || zm != mz) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Unexpected decomposition"); 1017*c4762a1bSJed Brown ierr = DMDAVecGetArray(da3,X3,&x);CHKERRQ(ierr); 1018*c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 1019*c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 1020*c4762a1bSJed Brown PetscReal u = PetscRealPart(x[i][j][zm-1].u); 1021*c4762a1bSJed Brown RangeUpdate(&umin,&umax,u); 1022*c4762a1bSJed Brown usum += u; 1023*c4762a1bSJed Brown } 1024*c4762a1bSJed Brown } 1025*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da3,X3,&x);CHKERRQ(ierr); 1026*c4762a1bSJed Brown ierr = DMCompositeRestoreAccess(pack,X,&X3,&X2);CHKERRQ(ierr); 1027*c4762a1bSJed Brown 1028*c4762a1bSJed Brown ierr = MPI_Allreduce(&umin,min,1,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)da3));CHKERRQ(ierr); 1029*c4762a1bSJed Brown ierr = MPI_Allreduce(&umax,max,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)da3));CHKERRQ(ierr); 1030*c4762a1bSJed Brown ierr = MPI_Allreduce(&usum,&gusum,1,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)da3));CHKERRQ(ierr); 1031*c4762a1bSJed Brown *mean = PetscRealPart(gusum) / (mx*my); 1032*c4762a1bSJed Brown PetscFunctionReturn(0); 1033*c4762a1bSJed Brown } 1034*c4762a1bSJed Brown 1035*c4762a1bSJed Brown static PetscErrorCode THISolveStatistics(THI thi,TS ts,PetscInt coarsened,const char name[]) 1036*c4762a1bSJed Brown { 1037*c4762a1bSJed Brown MPI_Comm comm; 1038*c4762a1bSJed Brown DM pack; 1039*c4762a1bSJed Brown Vec X,X3,X2; 1040*c4762a1bSJed Brown PetscErrorCode ierr; 1041*c4762a1bSJed Brown 1042*c4762a1bSJed Brown PetscFunctionBeginUser; 1043*c4762a1bSJed Brown ierr = PetscObjectGetComm((PetscObject)thi,&comm);CHKERRQ(ierr); 1044*c4762a1bSJed Brown ierr = TSGetDM(ts,&pack);CHKERRQ(ierr); 1045*c4762a1bSJed Brown ierr = TSGetSolution(ts,&X);CHKERRQ(ierr); 1046*c4762a1bSJed Brown ierr = DMCompositeGetAccess(pack,X,&X3,&X2);CHKERRQ(ierr); 1047*c4762a1bSJed Brown ierr = PetscPrintf(comm,"Solution statistics after solve: %s\n",name);CHKERRQ(ierr); 1048*c4762a1bSJed Brown { 1049*c4762a1bSJed Brown PetscInt its,lits; 1050*c4762a1bSJed Brown SNESConvergedReason reason; 1051*c4762a1bSJed Brown SNES snes; 1052*c4762a1bSJed Brown ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 1053*c4762a1bSJed Brown ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr); 1054*c4762a1bSJed Brown ierr = SNESGetConvergedReason(snes,&reason);CHKERRQ(ierr); 1055*c4762a1bSJed Brown ierr = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr); 1056*c4762a1bSJed Brown ierr = PetscPrintf(comm,"%s: Number of SNES iterations = %d, total linear iterations = %d\n",SNESConvergedReasons[reason],its,lits);CHKERRQ(ierr); 1057*c4762a1bSJed Brown } 1058*c4762a1bSJed Brown { 1059*c4762a1bSJed Brown PetscReal nrm2,tmin[3]={1e100,1e100,1e100},tmax[3]={-1e100,-1e100,-1e100},min[3],max[3]; 1060*c4762a1bSJed Brown PetscInt i,j,m; 1061*c4762a1bSJed Brown PetscScalar *x; 1062*c4762a1bSJed Brown ierr = VecNorm(X3,NORM_2,&nrm2);CHKERRQ(ierr); 1063*c4762a1bSJed Brown ierr = VecGetLocalSize(X3,&m);CHKERRQ(ierr); 1064*c4762a1bSJed Brown ierr = VecGetArray(X3,&x);CHKERRQ(ierr); 1065*c4762a1bSJed Brown for (i=0; i<m; i+=2) { 1066*c4762a1bSJed Brown PetscReal u = PetscRealPart(x[i]),v = PetscRealPart(x[i+1]),c = PetscSqrtReal(u*u+v*v); 1067*c4762a1bSJed Brown tmin[0] = PetscMin(u,tmin[0]); 1068*c4762a1bSJed Brown tmin[1] = PetscMin(v,tmin[1]); 1069*c4762a1bSJed Brown tmin[2] = PetscMin(c,tmin[2]); 1070*c4762a1bSJed Brown tmax[0] = PetscMax(u,tmax[0]); 1071*c4762a1bSJed Brown tmax[1] = PetscMax(v,tmax[1]); 1072*c4762a1bSJed Brown tmax[2] = PetscMax(c,tmax[2]); 1073*c4762a1bSJed Brown } 1074*c4762a1bSJed Brown ierr = VecRestoreArray(X,&x);CHKERRQ(ierr); 1075*c4762a1bSJed Brown ierr = MPI_Allreduce(tmin,min,3,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)thi));CHKERRQ(ierr); 1076*c4762a1bSJed Brown ierr = MPI_Allreduce(tmax,max,3,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)thi));CHKERRQ(ierr); 1077*c4762a1bSJed Brown /* Dimensionalize to meters/year */ 1078*c4762a1bSJed Brown nrm2 *= thi->units->year / thi->units->meter; 1079*c4762a1bSJed Brown for (j=0; j<3; j++) { 1080*c4762a1bSJed Brown min[j] *= thi->units->year / thi->units->meter; 1081*c4762a1bSJed Brown max[j] *= thi->units->year / thi->units->meter; 1082*c4762a1bSJed Brown } 1083*c4762a1bSJed Brown ierr = PetscPrintf(comm,"|X|_2 %g u in [%g, %g] v in [%g, %g] c in [%g, %g] \n",nrm2,min[0],max[0],min[1],max[1],min[2],max[2]);CHKERRQ(ierr); 1084*c4762a1bSJed Brown { 1085*c4762a1bSJed Brown PetscReal umin,umax,umean; 1086*c4762a1bSJed Brown ierr = THISurfaceStatistics(pack,X,&umin,&umax,&umean);CHKERRQ(ierr); 1087*c4762a1bSJed Brown umin *= thi->units->year / thi->units->meter; 1088*c4762a1bSJed Brown umax *= thi->units->year / thi->units->meter; 1089*c4762a1bSJed Brown umean *= thi->units->year / thi->units->meter; 1090*c4762a1bSJed Brown ierr = PetscPrintf(comm,"Surface statistics: u in [%12.6e, %12.6e] mean %12.6e\n",umin,umax,umean);CHKERRQ(ierr); 1091*c4762a1bSJed Brown } 1092*c4762a1bSJed Brown /* These values stay nondimensional */ 1093*c4762a1bSJed Brown ierr = PetscPrintf(comm,"Global eta range [%g, %g], converged range [%g, %g]\n",thi->eta.min,thi->eta.max,thi->eta.cmin,thi->eta.cmax);CHKERRQ(ierr); 1094*c4762a1bSJed Brown ierr = PetscPrintf(comm,"Global beta2 range [%g, %g], converged range [%g, %g]\n",thi->beta2.min,thi->beta2.max,thi->beta2.cmin,thi->beta2.cmax);CHKERRQ(ierr); 1095*c4762a1bSJed Brown } 1096*c4762a1bSJed Brown ierr = PetscPrintf(comm,"\n");CHKERRQ(ierr); 1097*c4762a1bSJed Brown ierr = DMCompositeRestoreAccess(pack,X,&X3,&X2);CHKERRQ(ierr); 1098*c4762a1bSJed Brown PetscFunctionReturn(0); 1099*c4762a1bSJed Brown } 1100*c4762a1bSJed Brown 1101*c4762a1bSJed Brown static inline PetscInt DMDALocalIndex3D(DMDALocalInfo *info,PetscInt i,PetscInt j,PetscInt k) 1102*c4762a1bSJed Brown {return ((i-info->gzs)*info->gym + (j-info->gys))*info->gxm + (k-info->gxs);} 1103*c4762a1bSJed Brown static inline PetscInt DMDALocalIndex2D(DMDALocalInfo *info,PetscInt i,PetscInt j) 1104*c4762a1bSJed Brown {return (i-info->gzs)*info->gym + (j-info->gys);} 1105*c4762a1bSJed Brown 1106*c4762a1bSJed Brown static PetscErrorCode THIJacobianLocal_Momentum(DMDALocalInfo *info,const Node ***x,const PrmNode **prm,Mat B,Mat Bcpl,THI thi) 1107*c4762a1bSJed Brown { 1108*c4762a1bSJed Brown PetscInt xs,ys,xm,ym,zm,i,j,k,q,l,ll; 1109*c4762a1bSJed Brown PetscReal hx,hy; 1110*c4762a1bSJed Brown PetscErrorCode ierr; 1111*c4762a1bSJed Brown 1112*c4762a1bSJed Brown PetscFunctionBeginUser; 1113*c4762a1bSJed Brown xs = info->zs; 1114*c4762a1bSJed Brown ys = info->ys; 1115*c4762a1bSJed Brown xm = info->zm; 1116*c4762a1bSJed Brown ym = info->ym; 1117*c4762a1bSJed Brown zm = info->xm; 1118*c4762a1bSJed Brown hx = thi->Lx / info->mz; 1119*c4762a1bSJed Brown hy = thi->Ly / info->my; 1120*c4762a1bSJed Brown 1121*c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 1122*c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 1123*c4762a1bSJed Brown PrmNode pn[4],dpn[4][2]; 1124*c4762a1bSJed Brown QuadExtract(prm,i,j,pn); 1125*c4762a1bSJed Brown ierr = QuadComputeGrad4(QuadQDeriv,hx,hy,pn,dpn);CHKERRQ(ierr); 1126*c4762a1bSJed Brown for (k=0; k<zm-1; k++) { 1127*c4762a1bSJed Brown Node n[8]; 1128*c4762a1bSJed Brown PetscReal zn[8],etabase = 0; 1129*c4762a1bSJed Brown PetscScalar Ke[8*NODE_SIZE][8*NODE_SIZE],Kcpl[8*NODE_SIZE][4*PRMNODE_SIZE]; 1130*c4762a1bSJed Brown PetscInt ls = 0; 1131*c4762a1bSJed Brown 1132*c4762a1bSJed Brown PrmHexGetZ(pn,k,zm,zn); 1133*c4762a1bSJed Brown HexExtract(x,i,j,k,n); 1134*c4762a1bSJed Brown ierr = PetscMemzero(Ke,sizeof(Ke));CHKERRQ(ierr); 1135*c4762a1bSJed Brown ierr = PetscMemzero(Kcpl,sizeof(Kcpl));CHKERRQ(ierr); 1136*c4762a1bSJed Brown if (thi->no_slip && k == 0) { 1137*c4762a1bSJed Brown for (l=0; l<4; l++) n[l].u = n[l].v = 0; 1138*c4762a1bSJed Brown ls = 4; 1139*c4762a1bSJed Brown } 1140*c4762a1bSJed Brown for (q=0; q<8; q++) { 1141*c4762a1bSJed Brown PetscReal dz[3],phi[8],dphi[8][3],jw,eta,deta; 1142*c4762a1bSJed Brown PetscScalar du[3],dv[3],u,v; 1143*c4762a1bSJed Brown HexGrad(HexQDeriv[q],zn,dz); 1144*c4762a1bSJed Brown HexComputeGeometry(q,hx,hy,dz,phi,dphi,&jw); 1145*c4762a1bSJed Brown PointwiseNonlinearity(thi,n,phi,dphi,&u,&v,du,dv,&eta,&deta); 1146*c4762a1bSJed Brown jw /= thi->rhog; /* residuals are scaled by this factor */ 1147*c4762a1bSJed Brown if (q == 0) etabase = eta; 1148*c4762a1bSJed Brown for (l=ls; l<8; l++) { /* test functions */ 1149*c4762a1bSJed Brown const PetscReal pp=phi[l],*restrict dp = dphi[l]; 1150*c4762a1bSJed Brown for (ll=ls; ll<8; ll++) { /* trial functions */ 1151*c4762a1bSJed Brown const PetscReal *restrict dpl = dphi[ll]; 1152*c4762a1bSJed Brown PetscScalar dgdu,dgdv; 1153*c4762a1bSJed Brown dgdu = 2.*du[0]*dpl[0] + dv[1]*dpl[0] + 0.5*(du[1]+dv[0])*dpl[1] + 0.5*du[2]*dpl[2]; 1154*c4762a1bSJed Brown dgdv = 2.*dv[1]*dpl[1] + du[0]*dpl[1] + 0.5*(du[1]+dv[0])*dpl[0] + 0.5*dv[2]*dpl[2]; 1155*c4762a1bSJed Brown /* Picard part */ 1156*c4762a1bSJed Brown 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]; 1157*c4762a1bSJed Brown Ke[l*2+0][ll*2+1] += dp[0]*jw*eta*2.*dpl[1] + dp[1]*jw*eta*dpl[0]; 1158*c4762a1bSJed Brown Ke[l*2+1][ll*2+0] += dp[1]*jw*eta*2.*dpl[0] + dp[0]*jw*eta*dpl[1]; 1159*c4762a1bSJed Brown 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]; 1160*c4762a1bSJed Brown /* extra Newton terms */ 1161*c4762a1bSJed Brown 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]; 1162*c4762a1bSJed Brown 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]; 1163*c4762a1bSJed Brown 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]; 1164*c4762a1bSJed Brown 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]; 1165*c4762a1bSJed Brown /* inertial part */ 1166*c4762a1bSJed Brown Ke[l*2+0][ll*2+0] += pp*jw*thi->inertia*pp; 1167*c4762a1bSJed Brown Ke[l*2+1][ll*2+1] += pp*jw*thi->inertia*pp; 1168*c4762a1bSJed Brown } 1169*c4762a1bSJed Brown for (ll=0; ll<4; ll++) { /* Trial functions for surface/bed */ 1170*c4762a1bSJed Brown const PetscReal dpl[] = {QuadQDeriv[q%4][ll][0]/hx, QuadQDeriv[q%4][ll][1]/hy}; /* surface = h + b */ 1171*c4762a1bSJed Brown Kcpl[FieldIndex(Node,l,u)][FieldIndex(PrmNode,ll,h)] += pp*jw*thi->rhog*dpl[0]; 1172*c4762a1bSJed Brown Kcpl[FieldIndex(Node,l,u)][FieldIndex(PrmNode,ll,b)] += pp*jw*thi->rhog*dpl[0]; 1173*c4762a1bSJed Brown Kcpl[FieldIndex(Node,l,v)][FieldIndex(PrmNode,ll,h)] += pp*jw*thi->rhog*dpl[1]; 1174*c4762a1bSJed Brown Kcpl[FieldIndex(Node,l,v)][FieldIndex(PrmNode,ll,b)] += pp*jw*thi->rhog*dpl[1]; 1175*c4762a1bSJed Brown } 1176*c4762a1bSJed Brown } 1177*c4762a1bSJed Brown } 1178*c4762a1bSJed Brown if (k == 0) { /* on a bottom face */ 1179*c4762a1bSJed Brown if (thi->no_slip) { 1180*c4762a1bSJed Brown const PetscReal hz = PetscRealPart(pn[0].h)/(zm-1); 1181*c4762a1bSJed Brown 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); 1182*c4762a1bSJed Brown Ke[0][0] = thi->dirichlet_scale*diagu; 1183*c4762a1bSJed Brown Ke[0][1] = 0; 1184*c4762a1bSJed Brown Ke[1][0] = 0; 1185*c4762a1bSJed Brown Ke[1][1] = thi->dirichlet_scale*diagv; 1186*c4762a1bSJed Brown } else { 1187*c4762a1bSJed Brown for (q=0; q<4; q++) { /* We remove the explicit scaling by 1/rhog because beta2 already has that scaling to be O(1) */ 1188*c4762a1bSJed Brown const PetscReal jw = 0.25*hx*hy,*phi = QuadQInterp[q]; 1189*c4762a1bSJed Brown PetscScalar u =0,v=0,rbeta2=0; 1190*c4762a1bSJed Brown PetscReal beta2,dbeta2; 1191*c4762a1bSJed Brown for (l=0; l<4; l++) { 1192*c4762a1bSJed Brown u += phi[l]*n[l].u; 1193*c4762a1bSJed Brown v += phi[l]*n[l].v; 1194*c4762a1bSJed Brown rbeta2 += phi[l]*pn[l].beta2; 1195*c4762a1bSJed Brown } 1196*c4762a1bSJed Brown THIFriction(thi,PetscRealPart(rbeta2),PetscRealPart(u*u+v*v)/2,&beta2,&dbeta2); 1197*c4762a1bSJed Brown for (l=0; l<4; l++) { 1198*c4762a1bSJed Brown const PetscReal pp = phi[l]; 1199*c4762a1bSJed Brown for (ll=0; ll<4; ll++) { 1200*c4762a1bSJed Brown const PetscReal ppl = phi[ll]; 1201*c4762a1bSJed Brown Ke[l*2+0][ll*2+0] += pp*jw*beta2*ppl + pp*jw*dbeta2*u*u*ppl; 1202*c4762a1bSJed Brown Ke[l*2+0][ll*2+1] += pp*jw*dbeta2*u*v*ppl; 1203*c4762a1bSJed Brown Ke[l*2+1][ll*2+0] += pp*jw*dbeta2*v*u*ppl; 1204*c4762a1bSJed Brown Ke[l*2+1][ll*2+1] += pp*jw*beta2*ppl + pp*jw*dbeta2*v*v*ppl; 1205*c4762a1bSJed Brown } 1206*c4762a1bSJed Brown } 1207*c4762a1bSJed Brown } 1208*c4762a1bSJed Brown } 1209*c4762a1bSJed Brown } 1210*c4762a1bSJed Brown { 1211*c4762a1bSJed Brown const PetscInt rc3blocked[8] = { 1212*c4762a1bSJed Brown DMDALocalIndex3D(info,i+0,j+0,k+0), 1213*c4762a1bSJed Brown DMDALocalIndex3D(info,i+1,j+0,k+0), 1214*c4762a1bSJed Brown DMDALocalIndex3D(info,i+1,j+1,k+0), 1215*c4762a1bSJed Brown DMDALocalIndex3D(info,i+0,j+1,k+0), 1216*c4762a1bSJed Brown DMDALocalIndex3D(info,i+0,j+0,k+1), 1217*c4762a1bSJed Brown DMDALocalIndex3D(info,i+1,j+0,k+1), 1218*c4762a1bSJed Brown DMDALocalIndex3D(info,i+1,j+1,k+1), 1219*c4762a1bSJed Brown DMDALocalIndex3D(info,i+0,j+1,k+1) 1220*c4762a1bSJed Brown },col2blocked[PRMNODE_SIZE*4] = { 1221*c4762a1bSJed Brown DMDALocalIndex2D(info,i+0,j+0), 1222*c4762a1bSJed Brown DMDALocalIndex2D(info,i+1,j+0), 1223*c4762a1bSJed Brown DMDALocalIndex2D(info,i+1,j+1), 1224*c4762a1bSJed Brown DMDALocalIndex2D(info,i+0,j+1) 1225*c4762a1bSJed Brown }; 1226*c4762a1bSJed Brown #if !defined COMPUTE_LOWER_TRIANGULAR /* fill in lower-triangular part, this is really cheap compared to computing the entries */ 1227*c4762a1bSJed Brown for (l=0; l<8; l++) { 1228*c4762a1bSJed Brown for (ll=l+1; ll<8; ll++) { 1229*c4762a1bSJed Brown Ke[ll*2+0][l*2+0] = Ke[l*2+0][ll*2+0]; 1230*c4762a1bSJed Brown Ke[ll*2+1][l*2+0] = Ke[l*2+0][ll*2+1]; 1231*c4762a1bSJed Brown Ke[ll*2+0][l*2+1] = Ke[l*2+1][ll*2+0]; 1232*c4762a1bSJed Brown Ke[ll*2+1][l*2+1] = Ke[l*2+1][ll*2+1]; 1233*c4762a1bSJed Brown } 1234*c4762a1bSJed Brown } 1235*c4762a1bSJed Brown #endif 1236*c4762a1bSJed Brown ierr = MatSetValuesBlockedLocal(B,8,rc3blocked,8,rc3blocked,&Ke[0][0],ADD_VALUES);CHKERRQ(ierr); /* velocity-velocity coupling can use blocked insertion */ 1237*c4762a1bSJed Brown { /* The off-diagonal part cannot (yet) */ 1238*c4762a1bSJed Brown PetscInt row3scalar[NODE_SIZE*8],col2scalar[PRMNODE_SIZE*4]; 1239*c4762a1bSJed Brown for (l=0; l<8; l++) for (ll=0; ll<NODE_SIZE; ll++) row3scalar[l*NODE_SIZE+ll] = rc3blocked[l]*NODE_SIZE+ll; 1240*c4762a1bSJed Brown for (l=0; l<4; l++) for (ll=0; ll<PRMNODE_SIZE; ll++) col2scalar[l*PRMNODE_SIZE+ll] = col2blocked[l]*PRMNODE_SIZE+ll; 1241*c4762a1bSJed Brown ierr = MatSetValuesLocal(Bcpl,8*NODE_SIZE,row3scalar,4*PRMNODE_SIZE,col2scalar,&Kcpl[0][0],ADD_VALUES);CHKERRQ(ierr); 1242*c4762a1bSJed Brown } 1243*c4762a1bSJed Brown } 1244*c4762a1bSJed Brown } 1245*c4762a1bSJed Brown } 1246*c4762a1bSJed Brown } 1247*c4762a1bSJed Brown PetscFunctionReturn(0); 1248*c4762a1bSJed Brown } 1249*c4762a1bSJed Brown 1250*c4762a1bSJed Brown static PetscErrorCode THIJacobianLocal_2D(DMDALocalInfo *info,const Node ***x3,const PrmNode **x2,const PrmNode **xdot2,PetscReal a,Mat B22,Mat B21,THI thi) 1251*c4762a1bSJed Brown { 1252*c4762a1bSJed Brown PetscErrorCode ierr; 1253*c4762a1bSJed Brown PetscInt xs,ys,xm,ym,zm,i,j,k; 1254*c4762a1bSJed Brown 1255*c4762a1bSJed Brown PetscFunctionBeginUser; 1256*c4762a1bSJed Brown xs = info->zs; 1257*c4762a1bSJed Brown ys = info->ys; 1258*c4762a1bSJed Brown xm = info->zm; 1259*c4762a1bSJed Brown ym = info->ym; 1260*c4762a1bSJed Brown zm = info->xm; 1261*c4762a1bSJed Brown 1262*c4762a1bSJed Brown if (zm > 1024) SETERRQ(((PetscObject)info->da)->comm,PETSC_ERR_SUP,"Need to allocate more space"); 1263*c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 1264*c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 1265*c4762a1bSJed Brown { /* Self-coupling */ 1266*c4762a1bSJed Brown const PetscInt row[] = {DMDALocalIndex2D(info,i,j)}; 1267*c4762a1bSJed Brown const PetscInt col[] = {DMDALocalIndex2D(info,i,j)}; 1268*c4762a1bSJed Brown const PetscScalar vals[] = { 1269*c4762a1bSJed Brown a,0,0, 1270*c4762a1bSJed Brown 0,a,0, 1271*c4762a1bSJed Brown 0,0,a 1272*c4762a1bSJed Brown }; 1273*c4762a1bSJed Brown ierr = MatSetValuesBlockedLocal(B22,1,row,1,col,vals,INSERT_VALUES);CHKERRQ(ierr); 1274*c4762a1bSJed Brown } 1275*c4762a1bSJed Brown for (k=0; k<zm; k++) { /* Coupling to velocity problem */ 1276*c4762a1bSJed Brown /* Use a cheaper quadrature than for residual evaluation, because it is much sparser */ 1277*c4762a1bSJed Brown const PetscInt row[] = {FieldIndex(PrmNode,DMDALocalIndex2D(info,i,j),h)}; 1278*c4762a1bSJed Brown const PetscInt cols[] = { 1279*c4762a1bSJed Brown FieldIndex(Node,DMDALocalIndex3D(info,i-1,j,k),u), 1280*c4762a1bSJed Brown FieldIndex(Node,DMDALocalIndex3D(info,i ,j,k),u), 1281*c4762a1bSJed Brown FieldIndex(Node,DMDALocalIndex3D(info,i+1,j,k),u), 1282*c4762a1bSJed Brown FieldIndex(Node,DMDALocalIndex3D(info,i,j-1,k),v), 1283*c4762a1bSJed Brown FieldIndex(Node,DMDALocalIndex3D(info,i,j ,k),v), 1284*c4762a1bSJed Brown FieldIndex(Node,DMDALocalIndex3D(info,i,j+1,k),v) 1285*c4762a1bSJed Brown }; 1286*c4762a1bSJed Brown const PetscScalar 1287*c4762a1bSJed Brown w = (k && k<zm-1) ? 0.5 : 0.25, 1288*c4762a1bSJed Brown hW = w*(x2[i-1][j ].h+x2[i ][j ].h)/(zm-1.), 1289*c4762a1bSJed Brown hE = w*(x2[i ][j ].h+x2[i+1][j ].h)/(zm-1.), 1290*c4762a1bSJed Brown hS = w*(x2[i ][j-1].h+x2[i ][j ].h)/(zm-1.), 1291*c4762a1bSJed Brown hN = w*(x2[i ][j ].h+x2[i ][j+1].h)/(zm-1.); 1292*c4762a1bSJed Brown PetscScalar *vals, 1293*c4762a1bSJed Brown vals_upwind[] = {((PetscRealPart(x3[i][j][k].u) > 0) ? -hW : 0), 1294*c4762a1bSJed Brown ((PetscRealPart(x3[i][j][k].u) > 0) ? +hE : -hW), 1295*c4762a1bSJed Brown ((PetscRealPart(x3[i][j][k].u) > 0) ? 0 : +hE), 1296*c4762a1bSJed Brown ((PetscRealPart(x3[i][j][k].v) > 0) ? -hS : 0), 1297*c4762a1bSJed Brown ((PetscRealPart(x3[i][j][k].v) > 0) ? +hN : -hS), 1298*c4762a1bSJed Brown ((PetscRealPart(x3[i][j][k].v) > 0) ? 0 : +hN)}, 1299*c4762a1bSJed Brown vals_centered[] = {-0.5*hW, 0.5*(-hW+hE), 0.5*hE, 1300*c4762a1bSJed Brown -0.5*hS, 0.5*(-hS+hN), 0.5*hN}; 1301*c4762a1bSJed Brown vals = 1 ? vals_upwind : vals_centered; 1302*c4762a1bSJed Brown if (k == 0) { 1303*c4762a1bSJed Brown Node derate; 1304*c4762a1bSJed Brown THIErosion(thi,&x3[i][j][0],NULL,&derate); 1305*c4762a1bSJed Brown vals[1] -= derate.u; 1306*c4762a1bSJed Brown vals[4] -= derate.v; 1307*c4762a1bSJed Brown } 1308*c4762a1bSJed Brown ierr = MatSetValuesLocal(B21,1,row,6,cols,vals,INSERT_VALUES);CHKERRQ(ierr); 1309*c4762a1bSJed Brown } 1310*c4762a1bSJed Brown } 1311*c4762a1bSJed Brown } 1312*c4762a1bSJed Brown PetscFunctionReturn(0); 1313*c4762a1bSJed Brown } 1314*c4762a1bSJed Brown 1315*c4762a1bSJed Brown static PetscErrorCode THIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,void *ctx) 1316*c4762a1bSJed Brown { 1317*c4762a1bSJed Brown PetscErrorCode ierr; 1318*c4762a1bSJed Brown THI thi = (THI)ctx; 1319*c4762a1bSJed Brown DM pack,da3,da2; 1320*c4762a1bSJed Brown Vec X3,X2,Xdot2; 1321*c4762a1bSJed Brown Mat B11,B12,B21,B22; 1322*c4762a1bSJed Brown DMDALocalInfo info3; 1323*c4762a1bSJed Brown IS *isloc; 1324*c4762a1bSJed Brown const Node ***x3; 1325*c4762a1bSJed Brown const PrmNode **x2,**xdot2; 1326*c4762a1bSJed Brown 1327*c4762a1bSJed Brown PetscFunctionBeginUser; 1328*c4762a1bSJed Brown ierr = TSGetDM(ts,&pack);CHKERRQ(ierr); 1329*c4762a1bSJed Brown ierr = DMCompositeGetEntries(pack,&da3,&da2);CHKERRQ(ierr); 1330*c4762a1bSJed Brown ierr = DMDAGetLocalInfo(da3,&info3);CHKERRQ(ierr); 1331*c4762a1bSJed Brown ierr = DMCompositeGetLocalVectors(pack,&X3,&X2);CHKERRQ(ierr); 1332*c4762a1bSJed Brown ierr = DMCompositeGetLocalVectors(pack,NULL,&Xdot2);CHKERRQ(ierr); 1333*c4762a1bSJed Brown ierr = DMCompositeScatter(pack,X,X3,X2);CHKERRQ(ierr); 1334*c4762a1bSJed Brown ierr = THIFixGhosts(thi,da3,da2,X3,X2);CHKERRQ(ierr); 1335*c4762a1bSJed Brown ierr = DMCompositeScatter(pack,Xdot,NULL,Xdot2);CHKERRQ(ierr); 1336*c4762a1bSJed Brown 1337*c4762a1bSJed Brown ierr = MatZeroEntries(B);CHKERRQ(ierr); 1338*c4762a1bSJed Brown 1339*c4762a1bSJed Brown ierr = DMCompositeGetLocalISs(pack,&isloc);CHKERRQ(ierr); 1340*c4762a1bSJed Brown ierr = MatGetLocalSubMatrix(B,isloc[0],isloc[0],&B11);CHKERRQ(ierr); 1341*c4762a1bSJed Brown ierr = MatGetLocalSubMatrix(B,isloc[0],isloc[1],&B12);CHKERRQ(ierr); 1342*c4762a1bSJed Brown ierr = MatGetLocalSubMatrix(B,isloc[1],isloc[0],&B21);CHKERRQ(ierr); 1343*c4762a1bSJed Brown ierr = MatGetLocalSubMatrix(B,isloc[1],isloc[1],&B22);CHKERRQ(ierr); 1344*c4762a1bSJed Brown 1345*c4762a1bSJed Brown ierr = DMDAVecGetArray(da3,X3,&x3);CHKERRQ(ierr); 1346*c4762a1bSJed Brown ierr = DMDAVecGetArray(da2,X2,&x2);CHKERRQ(ierr); 1347*c4762a1bSJed Brown ierr = DMDAVecGetArray(da2,Xdot2,&xdot2);CHKERRQ(ierr); 1348*c4762a1bSJed Brown 1349*c4762a1bSJed Brown ierr = THIJacobianLocal_Momentum(&info3,x3,x2,B11,B12,thi);CHKERRQ(ierr); 1350*c4762a1bSJed Brown 1351*c4762a1bSJed Brown /* Need to switch from ADD_VALUES to INSERT_VALUES */ 1352*c4762a1bSJed Brown ierr = MatAssemblyBegin(B,MAT_FLUSH_ASSEMBLY);CHKERRQ(ierr); 1353*c4762a1bSJed Brown ierr = MatAssemblyEnd(B,MAT_FLUSH_ASSEMBLY);CHKERRQ(ierr); 1354*c4762a1bSJed Brown 1355*c4762a1bSJed Brown ierr = THIJacobianLocal_2D(&info3,x3,x2,xdot2,a,B22,B21,thi);CHKERRQ(ierr); 1356*c4762a1bSJed Brown 1357*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da3,X3,&x3);CHKERRQ(ierr); 1358*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da2,X2,&x2);CHKERRQ(ierr); 1359*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da2,Xdot2,&xdot2);CHKERRQ(ierr); 1360*c4762a1bSJed Brown 1361*c4762a1bSJed Brown ierr = MatRestoreLocalSubMatrix(B,isloc[0],isloc[0],&B11);CHKERRQ(ierr); 1362*c4762a1bSJed Brown ierr = MatRestoreLocalSubMatrix(B,isloc[0],isloc[1],&B12);CHKERRQ(ierr); 1363*c4762a1bSJed Brown ierr = MatRestoreLocalSubMatrix(B,isloc[1],isloc[0],&B21);CHKERRQ(ierr); 1364*c4762a1bSJed Brown ierr = MatRestoreLocalSubMatrix(B,isloc[1],isloc[1],&B22);CHKERRQ(ierr); 1365*c4762a1bSJed Brown ierr = ISDestroy(&isloc[0]);CHKERRQ(ierr); 1366*c4762a1bSJed Brown ierr = ISDestroy(&isloc[1]);CHKERRQ(ierr); 1367*c4762a1bSJed Brown ierr = PetscFree(isloc);CHKERRQ(ierr); 1368*c4762a1bSJed Brown 1369*c4762a1bSJed Brown ierr = DMCompositeRestoreLocalVectors(pack,&X3,&X2);CHKERRQ(ierr); 1370*c4762a1bSJed Brown ierr = DMCompositeRestoreLocalVectors(pack,0,&Xdot2);CHKERRQ(ierr); 1371*c4762a1bSJed Brown 1372*c4762a1bSJed Brown ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1373*c4762a1bSJed Brown ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1374*c4762a1bSJed Brown if (A != B) { 1375*c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1376*c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1377*c4762a1bSJed Brown } 1378*c4762a1bSJed Brown if (thi->verbose) {ierr = THIMatrixStatistics(thi,B,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);} 1379*c4762a1bSJed Brown PetscFunctionReturn(0); 1380*c4762a1bSJed Brown } 1381*c4762a1bSJed Brown 1382*c4762a1bSJed Brown /* VTK's XML formats are so brain-dead that they can't handle multiple grids in the same file. Since the communication 1383*c4762a1bSJed Brown * can be shared between the two grids, we write two files at once, one for velocity (living on a 3D grid defined by 1384*c4762a1bSJed Brown * h=thickness and b=bed) and another for all properties living on the 2D grid. 1385*c4762a1bSJed Brown */ 1386*c4762a1bSJed Brown static PetscErrorCode THIDAVecView_VTK_XML(THI thi,DM pack,Vec X,const char filename[],const char filename2[]) 1387*c4762a1bSJed Brown { 1388*c4762a1bSJed Brown const PetscInt dof = NODE_SIZE,dof2 = PRMNODE_SIZE; 1389*c4762a1bSJed Brown Units units = thi->units; 1390*c4762a1bSJed Brown MPI_Comm comm; 1391*c4762a1bSJed Brown PetscErrorCode ierr; 1392*c4762a1bSJed Brown PetscViewer viewer3,viewer2; 1393*c4762a1bSJed Brown PetscMPIInt rank,size,tag,nn,nmax,nn2,nmax2; 1394*c4762a1bSJed Brown PetscInt mx,my,mz,r,range[6]; 1395*c4762a1bSJed Brown PetscScalar *x,*x2; 1396*c4762a1bSJed Brown DM da3,da2; 1397*c4762a1bSJed Brown Vec X3,X2; 1398*c4762a1bSJed Brown 1399*c4762a1bSJed Brown PetscFunctionBeginUser; 1400*c4762a1bSJed Brown ierr = PetscObjectGetComm((PetscObject)thi,&comm);CHKERRQ(ierr); 1401*c4762a1bSJed Brown ierr = DMCompositeGetEntries(pack,&da3,&da2);CHKERRQ(ierr); 1402*c4762a1bSJed Brown ierr = DMCompositeGetAccess(pack,X,&X3,&X2);CHKERRQ(ierr); 1403*c4762a1bSJed Brown ierr = DMDAGetInfo(da3,0, &mz,&my,&mx, 0,0,0, 0,0,0,0,0,0);CHKERRQ(ierr); 1404*c4762a1bSJed Brown ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1405*c4762a1bSJed Brown ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1406*c4762a1bSJed Brown ierr = PetscViewerASCIIOpen(comm,filename,&viewer3);CHKERRQ(ierr); 1407*c4762a1bSJed Brown ierr = PetscViewerASCIIOpen(comm,filename2,&viewer2);CHKERRQ(ierr); 1408*c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer3,"<VTKFile type=\"StructuredGrid\" version=\"0.1\" byte_order=\"LittleEndian\">\n");CHKERRQ(ierr); 1409*c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer2,"<VTKFile type=\"StructuredGrid\" version=\"0.1\" byte_order=\"LittleEndian\">\n");CHKERRQ(ierr); 1410*c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer3," <StructuredGrid WholeExtent=\"%d %d %d %d %d %d\">\n",0,mz-1,0,my-1,0,mx-1);CHKERRQ(ierr); 1411*c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer2," <StructuredGrid WholeExtent=\"%d %d %d %d %d %d\">\n",0,0,0,my-1,0,mx-1);CHKERRQ(ierr); 1412*c4762a1bSJed Brown 1413*c4762a1bSJed Brown ierr = DMDAGetCorners(da3,range,range+1,range+2,range+3,range+4,range+5);CHKERRQ(ierr); 1414*c4762a1bSJed Brown ierr = PetscMPIIntCast(range[3]*range[4]*range[5]*dof,&nn);CHKERRQ(ierr); 1415*c4762a1bSJed Brown ierr = MPI_Reduce(&nn,&nmax,1,MPI_INT,MPI_MAX,0,comm);CHKERRQ(ierr); 1416*c4762a1bSJed Brown ierr = PetscMPIIntCast(range[4]*range[5]*dof2,&nn2);CHKERRQ(ierr); 1417*c4762a1bSJed Brown ierr = MPI_Reduce(&nn2,&nmax2,1,MPI_INT,MPI_MAX,0,comm);CHKERRQ(ierr); 1418*c4762a1bSJed Brown tag = ((PetscObject)viewer3)->tag; 1419*c4762a1bSJed Brown ierr = VecGetArrayRead(X3,(const PetscScalar**)&x);CHKERRQ(ierr); 1420*c4762a1bSJed Brown ierr = VecGetArrayRead(X2,(const PetscScalar**)&x2);CHKERRQ(ierr); 1421*c4762a1bSJed Brown if (!rank) { 1422*c4762a1bSJed Brown PetscScalar *array,*array2; 1423*c4762a1bSJed Brown ierr = PetscMalloc2(nmax,&array,nmax2,&array2);CHKERRQ(ierr); 1424*c4762a1bSJed Brown for (r=0; r<size; r++) { 1425*c4762a1bSJed Brown PetscInt i,j,k,f,xs,xm,ys,ym,zs,zm; 1426*c4762a1bSJed Brown Node *y3; 1427*c4762a1bSJed Brown PetscScalar (*y2)[PRMNODE_SIZE]; 1428*c4762a1bSJed Brown MPI_Status status; 1429*c4762a1bSJed Brown 1430*c4762a1bSJed Brown if (r) { 1431*c4762a1bSJed Brown ierr = MPI_Recv(range,6,MPIU_INT,r,tag,comm,MPI_STATUS_IGNORE);CHKERRQ(ierr); 1432*c4762a1bSJed Brown } 1433*c4762a1bSJed Brown zs = range[0];ys = range[1];xs = range[2];zm = range[3];ym = range[4];xm = range[5]; 1434*c4762a1bSJed Brown if (xm*ym*zm*dof > nmax) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"should not happen"); 1435*c4762a1bSJed Brown if (r) { 1436*c4762a1bSJed Brown ierr = MPI_Recv(array,nmax,MPIU_SCALAR,r,tag,comm,&status);CHKERRQ(ierr); 1437*c4762a1bSJed Brown ierr = MPI_Get_count(&status,MPIU_SCALAR,&nn);CHKERRQ(ierr); 1438*c4762a1bSJed Brown if (nn != xm*ym*zm*dof) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"corrupt da3 send"); 1439*c4762a1bSJed Brown y3 = (Node*)array; 1440*c4762a1bSJed Brown ierr = MPI_Recv(array2,nmax2,MPIU_SCALAR,r,tag,comm,&status);CHKERRQ(ierr); 1441*c4762a1bSJed Brown ierr = MPI_Get_count(&status,MPIU_SCALAR,&nn2);CHKERRQ(ierr); 1442*c4762a1bSJed Brown if (nn2 != xm*ym*dof2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"corrupt da2 send"); 1443*c4762a1bSJed Brown y2 = (PetscScalar(*)[PRMNODE_SIZE])array2; 1444*c4762a1bSJed Brown } else { 1445*c4762a1bSJed Brown y3 = (Node*)x; 1446*c4762a1bSJed Brown y2 = (PetscScalar(*)[PRMNODE_SIZE])x2; 1447*c4762a1bSJed Brown } 1448*c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer3," <Piece Extent=\"%D %D %D %D %D %D\">\n",zs,zs+zm-1,ys,ys+ym-1,xs,xs+xm-1);CHKERRQ(ierr); 1449*c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer2," <Piece Extent=\"%d %d %D %D %D %D\">\n",0,0,ys,ys+ym-1,xs,xs+xm-1);CHKERRQ(ierr); 1450*c4762a1bSJed Brown 1451*c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer3," <Points>\n");CHKERRQ(ierr); 1452*c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer2," <Points>\n");CHKERRQ(ierr); 1453*c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer3," <DataArray type=\"Float32\" NumberOfComponents=\"3\" format=\"ascii\">\n");CHKERRQ(ierr); 1454*c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer2," <DataArray type=\"Float32\" NumberOfComponents=\"3\" format=\"ascii\">\n");CHKERRQ(ierr); 1455*c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 1456*c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 1457*c4762a1bSJed Brown PetscReal 1458*c4762a1bSJed Brown xx = thi->Lx*i/mx, 1459*c4762a1bSJed Brown yy = thi->Ly*j/my, 1460*c4762a1bSJed Brown b = PetscRealPart(y2[i*ym+j][FieldOffset(PrmNode,b)]), 1461*c4762a1bSJed Brown h = PetscRealPart(y2[i*ym+j][FieldOffset(PrmNode,h)]); 1462*c4762a1bSJed Brown for (k=zs; k<zs+zm; k++) { 1463*c4762a1bSJed Brown PetscReal zz = b + h*k/(mz-1.); 1464*c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer3,"%f %f %f\n",xx,yy,zz);CHKERRQ(ierr); 1465*c4762a1bSJed Brown } 1466*c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer2,"%f %f %f\n",xx,yy,(double)0.0);CHKERRQ(ierr); 1467*c4762a1bSJed Brown } 1468*c4762a1bSJed Brown } 1469*c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer3," </DataArray>\n");CHKERRQ(ierr); 1470*c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer2," </DataArray>\n");CHKERRQ(ierr); 1471*c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer3," </Points>\n");CHKERRQ(ierr); 1472*c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer2," </Points>\n");CHKERRQ(ierr); 1473*c4762a1bSJed Brown 1474*c4762a1bSJed Brown { /* Velocity and rank (3D) */ 1475*c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer3," <PointData>\n");CHKERRQ(ierr); 1476*c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer3," <DataArray type=\"Float32\" Name=\"velocity\" NumberOfComponents=\"3\" format=\"ascii\">\n");CHKERRQ(ierr); 1477*c4762a1bSJed Brown for (i=0; i<nn/dof; i++) { 1478*c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer3,"%f %f %f\n",PetscRealPart(y3[i].u)*units->year/units->meter,PetscRealPart(y3[i].v)*units->year/units->meter,0.0);CHKERRQ(ierr); 1479*c4762a1bSJed Brown } 1480*c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer3," </DataArray>\n");CHKERRQ(ierr); 1481*c4762a1bSJed Brown 1482*c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer3," <DataArray type=\"Int32\" Name=\"rank\" NumberOfComponents=\"1\" format=\"ascii\">\n");CHKERRQ(ierr); 1483*c4762a1bSJed Brown for (i=0; i<nn; i+=dof) { 1484*c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer3,"%D\n",r);CHKERRQ(ierr); 1485*c4762a1bSJed Brown } 1486*c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer3," </DataArray>\n");CHKERRQ(ierr); 1487*c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer3," </PointData>\n");CHKERRQ(ierr); 1488*c4762a1bSJed Brown } 1489*c4762a1bSJed Brown 1490*c4762a1bSJed Brown { /* 2D */ 1491*c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer2," <PointData>\n");CHKERRQ(ierr); 1492*c4762a1bSJed Brown for (f=0; f<PRMNODE_SIZE; f++) { 1493*c4762a1bSJed Brown const char *fieldname; 1494*c4762a1bSJed Brown ierr = DMDAGetFieldName(da2,f,&fieldname);CHKERRQ(ierr); 1495*c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer2," <DataArray type=\"Float32\" Name=\"%s\" format=\"ascii\">\n",fieldname);CHKERRQ(ierr); 1496*c4762a1bSJed Brown for (i=0; i<nn2/PRMNODE_SIZE; i++) { 1497*c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer2,"%g\n",y2[i][f]);CHKERRQ(ierr); 1498*c4762a1bSJed Brown } 1499*c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer2," </DataArray>\n");CHKERRQ(ierr); 1500*c4762a1bSJed Brown } 1501*c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer2," </PointData>\n");CHKERRQ(ierr); 1502*c4762a1bSJed Brown } 1503*c4762a1bSJed Brown 1504*c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer3," </Piece>\n");CHKERRQ(ierr); 1505*c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer2," </Piece>\n");CHKERRQ(ierr); 1506*c4762a1bSJed Brown } 1507*c4762a1bSJed Brown ierr = PetscFree2(array,array2);CHKERRQ(ierr); 1508*c4762a1bSJed Brown } else { 1509*c4762a1bSJed Brown ierr = MPI_Send(range,6,MPIU_INT,0,tag,comm);CHKERRQ(ierr); 1510*c4762a1bSJed Brown ierr = MPI_Send(x,nn,MPIU_SCALAR,0,tag,comm);CHKERRQ(ierr); 1511*c4762a1bSJed Brown ierr = MPI_Send(x2,nn2,MPIU_SCALAR,0,tag,comm);CHKERRQ(ierr); 1512*c4762a1bSJed Brown } 1513*c4762a1bSJed Brown ierr = VecRestoreArrayRead(X3,(const PetscScalar**)&x);CHKERRQ(ierr); 1514*c4762a1bSJed Brown ierr = VecRestoreArrayRead(X2,(const PetscScalar**)&x2);CHKERRQ(ierr); 1515*c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer3," </StructuredGrid>\n");CHKERRQ(ierr); 1516*c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer2," </StructuredGrid>\n");CHKERRQ(ierr); 1517*c4762a1bSJed Brown 1518*c4762a1bSJed Brown ierr = DMCompositeRestoreAccess(pack,X,&X3,&X2);CHKERRQ(ierr); 1519*c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer3,"</VTKFile>\n");CHKERRQ(ierr); 1520*c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer2,"</VTKFile>\n");CHKERRQ(ierr); 1521*c4762a1bSJed Brown ierr = PetscViewerDestroy(&viewer3);CHKERRQ(ierr); 1522*c4762a1bSJed Brown ierr = PetscViewerDestroy(&viewer2);CHKERRQ(ierr); 1523*c4762a1bSJed Brown PetscFunctionReturn(0); 1524*c4762a1bSJed Brown } 1525*c4762a1bSJed Brown 1526*c4762a1bSJed Brown static PetscErrorCode THITSMonitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx) 1527*c4762a1bSJed Brown { 1528*c4762a1bSJed Brown PetscErrorCode ierr; 1529*c4762a1bSJed Brown THI thi = (THI)ctx; 1530*c4762a1bSJed Brown DM pack; 1531*c4762a1bSJed Brown char filename3[PETSC_MAX_PATH_LEN],filename2[PETSC_MAX_PATH_LEN]; 1532*c4762a1bSJed Brown 1533*c4762a1bSJed Brown PetscFunctionBeginUser; 1534*c4762a1bSJed Brown if (step < 0) PetscFunctionReturn(0); /* negative one is used to indicate an interpolated solution */ 1535*c4762a1bSJed Brown ierr = PetscPrintf(PetscObjectComm((PetscObject)ts),"%3D: t=%g\n",step,(double)t);CHKERRQ(ierr); 1536*c4762a1bSJed Brown if (thi->monitor_interval && step % thi->monitor_interval) PetscFunctionReturn(0); 1537*c4762a1bSJed Brown ierr = TSGetDM(ts,&pack);CHKERRQ(ierr); 1538*c4762a1bSJed Brown ierr = PetscSNPrintf(filename3,sizeof(filename3),"%s-3d-%03d.vts",thi->monitor_basename,step);CHKERRQ(ierr); 1539*c4762a1bSJed Brown ierr = PetscSNPrintf(filename2,sizeof(filename2),"%s-2d-%03d.vts",thi->monitor_basename,step);CHKERRQ(ierr); 1540*c4762a1bSJed Brown ierr = THIDAVecView_VTK_XML(thi,pack,X,filename3,filename2);CHKERRQ(ierr); 1541*c4762a1bSJed Brown PetscFunctionReturn(0); 1542*c4762a1bSJed Brown } 1543*c4762a1bSJed Brown 1544*c4762a1bSJed Brown 1545*c4762a1bSJed Brown static PetscErrorCode THICreateDM3d(THI thi,DM *dm3d) 1546*c4762a1bSJed Brown { 1547*c4762a1bSJed Brown MPI_Comm comm; 1548*c4762a1bSJed Brown PetscInt M = 3,N = 3,P = 2; 1549*c4762a1bSJed Brown DM da; 1550*c4762a1bSJed Brown PetscErrorCode ierr; 1551*c4762a1bSJed Brown 1552*c4762a1bSJed Brown PetscFunctionBeginUser; 1553*c4762a1bSJed Brown ierr = PetscObjectGetComm((PetscObject)thi,&comm);CHKERRQ(ierr); 1554*c4762a1bSJed Brown ierr = PetscOptionsBegin(comm,NULL,"Grid resolution options","");CHKERRQ(ierr); 1555*c4762a1bSJed Brown { 1556*c4762a1bSJed Brown ierr = PetscOptionsInt("-M","Number of elements in x-direction on coarse level","",M,&M,NULL);CHKERRQ(ierr); 1557*c4762a1bSJed Brown N = M; 1558*c4762a1bSJed Brown ierr = PetscOptionsInt("-N","Number of elements in y-direction on coarse level (if different from M)","",N,&N,NULL);CHKERRQ(ierr); 1559*c4762a1bSJed Brown ierr = PetscOptionsInt("-P","Number of elements in z-direction on coarse level","",P,&P,NULL);CHKERRQ(ierr); 1560*c4762a1bSJed Brown } 1561*c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 1562*c4762a1bSJed Brown ierr = 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);CHKERRQ(ierr); 1563*c4762a1bSJed Brown ierr = DMSetFromOptions(da);CHKERRQ(ierr); 1564*c4762a1bSJed Brown ierr = DMSetUp(da);CHKERRQ(ierr); 1565*c4762a1bSJed Brown ierr = DMDASetFieldName(da,0,"x-velocity");CHKERRQ(ierr); 1566*c4762a1bSJed Brown ierr = DMDASetFieldName(da,1,"y-velocity");CHKERRQ(ierr); 1567*c4762a1bSJed Brown *dm3d = da; 1568*c4762a1bSJed Brown PetscFunctionReturn(0); 1569*c4762a1bSJed Brown } 1570*c4762a1bSJed Brown 1571*c4762a1bSJed Brown int main(int argc,char *argv[]) 1572*c4762a1bSJed Brown { 1573*c4762a1bSJed Brown MPI_Comm comm; 1574*c4762a1bSJed Brown DM pack,da3,da2; 1575*c4762a1bSJed Brown TS ts; 1576*c4762a1bSJed Brown THI thi; 1577*c4762a1bSJed Brown Vec X; 1578*c4762a1bSJed Brown Mat B; 1579*c4762a1bSJed Brown PetscInt i,steps; 1580*c4762a1bSJed Brown PetscReal ftime; 1581*c4762a1bSJed Brown PetscErrorCode ierr; 1582*c4762a1bSJed Brown 1583*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,0,help);if (ierr) return ierr; 1584*c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 1585*c4762a1bSJed Brown 1586*c4762a1bSJed Brown ierr = THICreate(comm,&thi);CHKERRQ(ierr); 1587*c4762a1bSJed Brown ierr = THICreateDM3d(thi,&da3);CHKERRQ(ierr); 1588*c4762a1bSJed Brown { 1589*c4762a1bSJed Brown PetscInt Mx,My,mx,my,s; 1590*c4762a1bSJed Brown DMDAStencilType st; 1591*c4762a1bSJed Brown ierr = DMDAGetInfo(da3,0, 0,&My,&Mx, 0,&my,&mx, 0,&s,0,0,0,&st);CHKERRQ(ierr); 1592*c4762a1bSJed Brown ierr = DMDACreate2d(PetscObjectComm((PetscObject)thi),DM_BOUNDARY_PERIODIC,DM_BOUNDARY_PERIODIC,st,My,Mx,my,mx,sizeof(PrmNode)/sizeof(PetscScalar),s,0,0,&da2);CHKERRQ(ierr); 1593*c4762a1bSJed Brown ierr = DMSetUp(da2);CHKERRQ(ierr); 1594*c4762a1bSJed Brown } 1595*c4762a1bSJed Brown 1596*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject)da3,"3D_Velocity");CHKERRQ(ierr); 1597*c4762a1bSJed Brown ierr = DMSetOptionsPrefix(da3,"f3d_");CHKERRQ(ierr); 1598*c4762a1bSJed Brown ierr = DMDASetFieldName(da3,0,"u");CHKERRQ(ierr); 1599*c4762a1bSJed Brown ierr = DMDASetFieldName(da3,1,"v");CHKERRQ(ierr); 1600*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject)da2,"2D_Fields");CHKERRQ(ierr); 1601*c4762a1bSJed Brown ierr = DMSetOptionsPrefix(da2,"f2d_");CHKERRQ(ierr); 1602*c4762a1bSJed Brown ierr = DMDASetFieldName(da2,0,"b");CHKERRQ(ierr); 1603*c4762a1bSJed Brown ierr = DMDASetFieldName(da2,1,"h");CHKERRQ(ierr); 1604*c4762a1bSJed Brown ierr = DMDASetFieldName(da2,2,"beta2");CHKERRQ(ierr); 1605*c4762a1bSJed Brown ierr = DMCompositeCreate(comm,&pack);CHKERRQ(ierr); 1606*c4762a1bSJed Brown ierr = DMCompositeAddDM(pack,da3);CHKERRQ(ierr); 1607*c4762a1bSJed Brown ierr = DMCompositeAddDM(pack,da2);CHKERRQ(ierr); 1608*c4762a1bSJed Brown ierr = DMDestroy(&da3);CHKERRQ(ierr); 1609*c4762a1bSJed Brown ierr = DMDestroy(&da2);CHKERRQ(ierr); 1610*c4762a1bSJed Brown ierr = DMSetUp(pack);CHKERRQ(ierr); 1611*c4762a1bSJed Brown ierr = DMCreateMatrix(pack,&B);CHKERRQ(ierr); 1612*c4762a1bSJed Brown ierr = MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr); 1613*c4762a1bSJed Brown ierr = MatSetOptionsPrefix(B,"thi_");CHKERRQ(ierr); 1614*c4762a1bSJed Brown 1615*c4762a1bSJed Brown for (i=0; i<thi->nlevels; i++) { 1616*c4762a1bSJed Brown PetscReal Lx = thi->Lx / thi->units->meter,Ly = thi->Ly / thi->units->meter,Lz = thi->Lz / thi->units->meter; 1617*c4762a1bSJed Brown PetscInt Mx,My,Mz; 1618*c4762a1bSJed Brown ierr = DMCompositeGetEntries(pack,&da3,&da2);CHKERRQ(ierr); 1619*c4762a1bSJed Brown ierr = DMDAGetInfo(da3,0, &Mz,&My,&Mx, 0,0,0, 0,0,0,0,0,0);CHKERRQ(ierr); 1620*c4762a1bSJed Brown ierr = PetscPrintf(PetscObjectComm((PetscObject)thi),"Level %D domain size (m) %8.2g x %8.2g x %8.2g, num elements %3d x %3d x %3d (%8d), size (m) %g x %g x %g\n",i,Lx,Ly,Lz,Mx,My,Mz,Mx*My*Mz,Lx/Mx,Ly/My,1000./(Mz-1));CHKERRQ(ierr); 1621*c4762a1bSJed Brown } 1622*c4762a1bSJed Brown 1623*c4762a1bSJed Brown ierr = DMCreateGlobalVector(pack,&X);CHKERRQ(ierr); 1624*c4762a1bSJed Brown ierr = THIInitial(thi,pack,X);CHKERRQ(ierr); 1625*c4762a1bSJed Brown 1626*c4762a1bSJed Brown ierr = TSCreate(comm,&ts);CHKERRQ(ierr); 1627*c4762a1bSJed Brown ierr = TSSetDM(ts,pack);CHKERRQ(ierr); 1628*c4762a1bSJed Brown ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr); 1629*c4762a1bSJed Brown ierr = TSMonitorSet(ts,THITSMonitor,thi,NULL);CHKERRQ(ierr); 1630*c4762a1bSJed Brown ierr = TSSetType(ts,TSTHETA);CHKERRQ(ierr); 1631*c4762a1bSJed Brown ierr = TSSetIFunction(ts,NULL,THIFunction,thi);CHKERRQ(ierr); 1632*c4762a1bSJed Brown ierr = TSSetIJacobian(ts,B,B,THIJacobian,thi);CHKERRQ(ierr); 1633*c4762a1bSJed Brown ierr = TSSetMaxTime(ts,10.0);CHKERRQ(ierr); 1634*c4762a1bSJed Brown ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr); 1635*c4762a1bSJed Brown ierr = TSSetSolution(ts,X);CHKERRQ(ierr); 1636*c4762a1bSJed Brown ierr = TSSetTimeStep(ts,1e-3);CHKERRQ(ierr); 1637*c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 1638*c4762a1bSJed Brown 1639*c4762a1bSJed Brown ierr = TSSolve(ts,X);CHKERRQ(ierr); 1640*c4762a1bSJed Brown ierr = TSGetSolveTime(ts,&ftime);CHKERRQ(ierr); 1641*c4762a1bSJed Brown ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr); 1642*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Steps %D final time %g\n",steps,(double)ftime);CHKERRQ(ierr); 1643*c4762a1bSJed Brown 1644*c4762a1bSJed Brown if (0) {ierr = THISolveStatistics(thi,ts,0,"Full");CHKERRQ(ierr);} 1645*c4762a1bSJed Brown 1646*c4762a1bSJed Brown { 1647*c4762a1bSJed Brown PetscBool flg; 1648*c4762a1bSJed Brown char filename[PETSC_MAX_PATH_LEN] = ""; 1649*c4762a1bSJed Brown ierr = PetscOptionsGetString(NULL,NULL,"-o",filename,sizeof(filename),&flg);CHKERRQ(ierr); 1650*c4762a1bSJed Brown if (flg) { 1651*c4762a1bSJed Brown ierr = THIDAVecView_VTK_XML(thi,pack,X,filename,NULL);CHKERRQ(ierr); 1652*c4762a1bSJed Brown } 1653*c4762a1bSJed Brown } 1654*c4762a1bSJed Brown 1655*c4762a1bSJed Brown ierr = VecDestroy(&X);CHKERRQ(ierr); 1656*c4762a1bSJed Brown ierr = MatDestroy(&B);CHKERRQ(ierr); 1657*c4762a1bSJed Brown ierr = DMDestroy(&pack);CHKERRQ(ierr); 1658*c4762a1bSJed Brown ierr = TSDestroy(&ts);CHKERRQ(ierr); 1659*c4762a1bSJed Brown ierr = THIDestroy(&thi);CHKERRQ(ierr); 1660*c4762a1bSJed Brown ierr = PetscFinalize(); 1661*c4762a1bSJed Brown return ierr; 1662*c4762a1bSJed Brown } 1663