1c4762a1bSJed Brown static const char help[] = "Toy hydrostatic ice flow with multigrid in 3D.\n\ 2c4762a1bSJed Brown \n\ 3c4762a1bSJed Brown Solves the hydrostatic (aka Blatter/Pattyn/First Order) equations for ice sheet flow\n\ 4c4762a1bSJed Brown using multigrid. The ice uses a power-law rheology with \"Glen\" exponent 3 (corresponds\n\ 5c4762a1bSJed Brown to p=4/3 in a p-Laplacian). The focus is on ISMIP-HOM experiments which assume periodic\n\ 6c4762a1bSJed Brown boundary conditions in the x- and y-directions.\n\ 7c4762a1bSJed Brown \n\ 8c4762a1bSJed Brown Equations are rescaled so that the domain size and solution are O(1), details of this scaling\n\ 9c4762a1bSJed Brown can be controlled by the options -units_meter, -units_second, and -units_kilogram.\n\ 10c4762a1bSJed Brown \n\ 11c4762a1bSJed Brown A VTK StructuredGrid output file can be written using the option -o filename.vts\n\ 12c4762a1bSJed Brown \n\n"; 13c4762a1bSJed Brown 14c4762a1bSJed Brown /* 15c4762a1bSJed Brown The equations for horizontal velocity (u,v) are 16c4762a1bSJed Brown 17c4762a1bSJed Brown - [eta (4 u_x + 2 v_y)]_x - [eta (u_y + v_x)]_y - [eta u_z]_z + rho g s_x = 0 18c4762a1bSJed Brown - [eta (4 v_y + 2 u_x)]_y - [eta (u_y + v_x)]_x - [eta v_z]_z + rho g s_y = 0 19c4762a1bSJed Brown 20c4762a1bSJed Brown where 21c4762a1bSJed Brown 22c4762a1bSJed Brown eta = B/2 (epsilon + gamma)^((p-2)/2) 23c4762a1bSJed Brown 24c4762a1bSJed Brown is the nonlinear effective viscosity with regularization epsilon and hardness parameter B, 25c4762a1bSJed Brown written in terms of the second invariant 26c4762a1bSJed Brown 27c4762a1bSJed 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 28c4762a1bSJed Brown 29c4762a1bSJed Brown The surface boundary conditions are the natural conditions. The basal boundary conditions 30c4762a1bSJed Brown are either no-slip, or Navier (linear) slip with spatially variant friction coefficient beta^2. 31c4762a1bSJed Brown 32c4762a1bSJed Brown In the code, the equations for (u,v) are multiplied through by 1/(rho g) so that residuals are O(1). 33c4762a1bSJed Brown 34c4762a1bSJed Brown The discretization is Q1 finite elements, managed by a DMDA. The grid is never distorted in the 35c4762a1bSJed Brown map (x,y) plane, but the bed and surface may be bumpy. This is handled as usual in FEM, through 36c4762a1bSJed Brown the Jacobian of the coordinate transformation from a reference element to the physical element. 37c4762a1bSJed Brown 38c4762a1bSJed Brown Since ice-flow is tightly coupled in the z-direction (within columns), the DMDA is managed 39c4762a1bSJed Brown specially so that columns are never distributed, and are always contiguous in memory. 40c4762a1bSJed Brown This amounts to reversing the meaning of X,Y,Z compared to the DMDA's internal interpretation, 41c4762a1bSJed Brown and then indexing as vec[i][j][k]. The exotic coarse spaces require 2D DMDAs which are made to 42c4762a1bSJed Brown use compatible domain decomposition relative to the 3D DMDAs. 43c4762a1bSJed Brown 44c4762a1bSJed Brown */ 45c4762a1bSJed Brown 46c4762a1bSJed Brown #include <petscts.h> 47c4762a1bSJed Brown #include <petscdm.h> 48c4762a1bSJed Brown #include <petscdmda.h> 49c4762a1bSJed Brown #include <petscdmcomposite.h> 50c4762a1bSJed Brown #include <ctype.h> /* toupper() */ 51c4762a1bSJed Brown #include <petsc/private/petscimpl.h> 52c4762a1bSJed Brown 53c4762a1bSJed Brown #if defined __SSE2__ 54c4762a1bSJed Brown # include <emmintrin.h> 55c4762a1bSJed Brown #endif 56c4762a1bSJed Brown 57c4762a1bSJed Brown /* The SSE2 kernels are only for PetscScalar=double on architectures that support it */ 58c4762a1bSJed Brown #define USE_SSE2_KERNELS (!defined NO_SSE2 \ 59c4762a1bSJed Brown && !defined PETSC_USE_COMPLEX \ 60c4762a1bSJed Brown && !defined PETSC_USE_REAL_SINGLE \ 61c4762a1bSJed Brown && defined __SSE2__) 62c4762a1bSJed Brown 63c4762a1bSJed Brown #if !defined __STDC_VERSION__ || __STDC_VERSION__ < 199901L 64c4762a1bSJed Brown # if defined __cplusplus /* C++ restrict is nonstandard and compilers have inconsistent rules about where it can be used */ 65c4762a1bSJed Brown # define restrict 66c4762a1bSJed Brown # else 67c4762a1bSJed Brown # define restrict PETSC_RESTRICT 68c4762a1bSJed Brown # endif 69c4762a1bSJed Brown #endif 70c4762a1bSJed Brown 71c4762a1bSJed Brown static PetscClassId THI_CLASSID; 72c4762a1bSJed Brown 73c4762a1bSJed Brown typedef enum {QUAD_GAUSS,QUAD_LOBATTO} QuadratureType; 74c4762a1bSJed Brown static const char *QuadratureTypes[] = {"gauss","lobatto","QuadratureType","QUAD_",0}; 75c4762a1bSJed Brown static const PetscReal HexQWeights[8] = {1,1,1,1,1,1,1,1}; 76c4762a1bSJed Brown static const PetscReal HexQNodes[] = {-0.57735026918962573, 0.57735026918962573}; 77c4762a1bSJed Brown #define G 0.57735026918962573 78c4762a1bSJed Brown #define H (0.5*(1.+G)) 79c4762a1bSJed Brown #define L (0.5*(1.-G)) 80c4762a1bSJed Brown #define M (-0.5) 81c4762a1bSJed Brown #define P (0.5) 82c4762a1bSJed Brown /* Special quadrature: Lobatto in horizontal, Gauss in vertical */ 83c4762a1bSJed Brown static const PetscReal HexQInterp_Lobatto[8][8] = {{H,0,0,0,L,0,0,0}, 84c4762a1bSJed Brown {0,H,0,0,0,L,0,0}, 85c4762a1bSJed Brown {0,0,H,0,0,0,L,0}, 86c4762a1bSJed Brown {0,0,0,H,0,0,0,L}, 87c4762a1bSJed Brown {L,0,0,0,H,0,0,0}, 88c4762a1bSJed Brown {0,L,0,0,0,H,0,0}, 89c4762a1bSJed Brown {0,0,L,0,0,0,H,0}, 90c4762a1bSJed Brown {0,0,0,L,0,0,0,H}}; 91c4762a1bSJed Brown static const PetscReal HexQDeriv_Lobatto[8][8][3] = { 92c4762a1bSJed 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} }, 93c4762a1bSJed 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} }, 94c4762a1bSJed 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} }, 95c4762a1bSJed 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}}, 96c4762a1bSJed 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} }, 97c4762a1bSJed 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} }, 98c4762a1bSJed 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} }, 99c4762a1bSJed 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}}}; 100c4762a1bSJed Brown /* Stanndard Gauss */ 101c4762a1bSJed 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}, 102c4762a1bSJed 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}, 103c4762a1bSJed 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}, 104c4762a1bSJed 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}, 105c4762a1bSJed 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}, 106c4762a1bSJed 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}, 107c4762a1bSJed 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}, 108c4762a1bSJed 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}}; 109c4762a1bSJed Brown static const PetscReal HexQDeriv_Gauss[8][8][3] = { 110c4762a1bSJed 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}}, 111c4762a1bSJed 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}}, 112c4762a1bSJed 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}}, 113c4762a1bSJed 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}}, 114c4762a1bSJed 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}}, 115c4762a1bSJed 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}}, 116c4762a1bSJed 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}}, 117c4762a1bSJed 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}}}; 118c4762a1bSJed Brown static const PetscReal (*HexQInterp)[8],(*HexQDeriv)[8][3]; 119c4762a1bSJed Brown /* Standard 2x2 Gauss quadrature for the bottom layer. */ 120c4762a1bSJed Brown static const PetscReal QuadQInterp[4][4] = {{H*H,L*H,L*L,H*L}, 121c4762a1bSJed Brown {L*H,H*H,H*L,L*L}, 122c4762a1bSJed Brown {L*L,H*L,H*H,L*H}, 123c4762a1bSJed Brown {H*L,L*L,L*H,H*H}}; 124c4762a1bSJed Brown static const PetscReal QuadQDeriv[4][4][2] = { 125c4762a1bSJed Brown {{M*H,M*H},{P*H,M*L},{P*L,P*L},{M*L,P*H}}, 126c4762a1bSJed Brown {{M*H,M*L},{P*H,M*H},{P*L,P*H},{M*L,P*L}}, 127c4762a1bSJed Brown {{M*L,M*L},{P*L,M*H},{P*H,P*H},{M*H,P*L}}, 128c4762a1bSJed Brown {{M*L,M*H},{P*L,M*L},{P*H,P*L},{M*H,P*H}}}; 129c4762a1bSJed Brown #undef G 130c4762a1bSJed Brown #undef H 131c4762a1bSJed Brown #undef L 132c4762a1bSJed Brown #undef M 133c4762a1bSJed Brown #undef P 134c4762a1bSJed Brown 135c4762a1bSJed Brown #define HexExtract(x,i,j,k,n) do { \ 136c4762a1bSJed Brown (n)[0] = (x)[i][j][k]; \ 137c4762a1bSJed Brown (n)[1] = (x)[i+1][j][k]; \ 138c4762a1bSJed Brown (n)[2] = (x)[i+1][j+1][k]; \ 139c4762a1bSJed Brown (n)[3] = (x)[i][j+1][k]; \ 140c4762a1bSJed Brown (n)[4] = (x)[i][j][k+1]; \ 141c4762a1bSJed Brown (n)[5] = (x)[i+1][j][k+1]; \ 142c4762a1bSJed Brown (n)[6] = (x)[i+1][j+1][k+1]; \ 143c4762a1bSJed Brown (n)[7] = (x)[i][j+1][k+1]; \ 144c4762a1bSJed Brown } while (0) 145c4762a1bSJed Brown 146c4762a1bSJed Brown #define HexExtractRef(x,i,j,k,n) do { \ 147c4762a1bSJed Brown (n)[0] = &(x)[i][j][k]; \ 148c4762a1bSJed Brown (n)[1] = &(x)[i+1][j][k]; \ 149c4762a1bSJed Brown (n)[2] = &(x)[i+1][j+1][k]; \ 150c4762a1bSJed Brown (n)[3] = &(x)[i][j+1][k]; \ 151c4762a1bSJed Brown (n)[4] = &(x)[i][j][k+1]; \ 152c4762a1bSJed Brown (n)[5] = &(x)[i+1][j][k+1]; \ 153c4762a1bSJed Brown (n)[6] = &(x)[i+1][j+1][k+1]; \ 154c4762a1bSJed Brown (n)[7] = &(x)[i][j+1][k+1]; \ 155c4762a1bSJed Brown } while (0) 156c4762a1bSJed Brown 157c4762a1bSJed Brown #define QuadExtract(x,i,j,n) do { \ 158c4762a1bSJed Brown (n)[0] = (x)[i][j]; \ 159c4762a1bSJed Brown (n)[1] = (x)[i+1][j]; \ 160c4762a1bSJed Brown (n)[2] = (x)[i+1][j+1]; \ 161c4762a1bSJed Brown (n)[3] = (x)[i][j+1]; \ 162c4762a1bSJed Brown } while (0) 163c4762a1bSJed Brown 164c4762a1bSJed Brown static PetscScalar Sqr(PetscScalar a) {return a*a;} 165c4762a1bSJed Brown 166c4762a1bSJed Brown static void HexGrad(const PetscReal dphi[][3],const PetscReal zn[],PetscReal dz[]) 167c4762a1bSJed Brown { 168c4762a1bSJed Brown PetscInt i; 169c4762a1bSJed Brown dz[0] = dz[1] = dz[2] = 0; 170c4762a1bSJed Brown for (i=0; i<8; i++) { 171c4762a1bSJed Brown dz[0] += dphi[i][0] * zn[i]; 172c4762a1bSJed Brown dz[1] += dphi[i][1] * zn[i]; 173c4762a1bSJed Brown dz[2] += dphi[i][2] * zn[i]; 174c4762a1bSJed Brown } 175c4762a1bSJed Brown } 176c4762a1bSJed Brown 177c4762a1bSJed Brown static void HexComputeGeometry(PetscInt q,PetscReal hx,PetscReal hy,const PetscReal dz[restrict],PetscReal phi[restrict],PetscReal dphi[restrict][3],PetscReal *restrict jw) 178c4762a1bSJed Brown { 179c4762a1bSJed Brown const PetscReal 180c4762a1bSJed Brown jac[3][3] = {{hx/2,0,0}, {0,hy/2,0}, {dz[0],dz[1],dz[2]}} 181c4762a1bSJed 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]}} 182c4762a1bSJed Brown ,jdet = jac[0][0]*jac[1][1]*jac[2][2]; 183c4762a1bSJed Brown PetscInt i; 184c4762a1bSJed Brown 185c4762a1bSJed Brown for (i=0; i<8; i++) { 186c4762a1bSJed Brown const PetscReal *dphir = HexQDeriv[q][i]; 187c4762a1bSJed Brown phi[i] = HexQInterp[q][i]; 188c4762a1bSJed Brown dphi[i][0] = dphir[0]*ijac[0][0] + dphir[1]*ijac[1][0] + dphir[2]*ijac[2][0]; 189c4762a1bSJed Brown dphi[i][1] = dphir[0]*ijac[0][1] + dphir[1]*ijac[1][1] + dphir[2]*ijac[2][1]; 190c4762a1bSJed Brown dphi[i][2] = dphir[0]*ijac[0][2] + dphir[1]*ijac[1][2] + dphir[2]*ijac[2][2]; 191c4762a1bSJed Brown } 192c4762a1bSJed Brown *jw = 1.0 * jdet; 193c4762a1bSJed Brown } 194c4762a1bSJed Brown 195c4762a1bSJed Brown typedef struct _p_THI *THI; 196c4762a1bSJed Brown typedef struct _n_Units *Units; 197c4762a1bSJed Brown 198c4762a1bSJed Brown typedef struct { 199c4762a1bSJed Brown PetscScalar u,v; 200c4762a1bSJed Brown } Node; 201c4762a1bSJed Brown 202c4762a1bSJed Brown typedef struct { 203c4762a1bSJed Brown PetscScalar b; /* bed */ 204c4762a1bSJed Brown PetscScalar h; /* thickness */ 205c4762a1bSJed Brown PetscScalar beta2; /* friction */ 206c4762a1bSJed Brown } PrmNode; 207c4762a1bSJed Brown 208c4762a1bSJed Brown #define FieldSize(ntype) ((PetscInt)(sizeof(ntype)/sizeof(PetscScalar))) 209c4762a1bSJed Brown #define FieldOffset(ntype,member) ((PetscInt)(offsetof(ntype,member)/sizeof(PetscScalar))) 210c4762a1bSJed Brown #define FieldIndex(ntype,i,member) ((PetscInt)((i)*FieldSize(ntype) + FieldOffset(ntype,member))) 211c4762a1bSJed Brown #define NODE_SIZE FieldSize(Node) 212c4762a1bSJed Brown #define PRMNODE_SIZE FieldSize(PrmNode) 213c4762a1bSJed Brown 214c4762a1bSJed Brown typedef struct { 215c4762a1bSJed Brown PetscReal min,max,cmin,cmax; 216c4762a1bSJed Brown } PRange; 217c4762a1bSJed Brown 218c4762a1bSJed Brown struct _p_THI { 219c4762a1bSJed Brown PETSCHEADER(int); 220c4762a1bSJed Brown void (*initialize)(THI,PetscReal x,PetscReal y,PrmNode *p); 221c4762a1bSJed Brown PetscInt nlevels; 222c4762a1bSJed Brown PetscInt zlevels; 223c4762a1bSJed Brown PetscReal Lx,Ly,Lz; /* Model domain */ 224c4762a1bSJed Brown PetscReal alpha; /* Bed angle */ 225c4762a1bSJed Brown Units units; 226c4762a1bSJed Brown PetscReal dirichlet_scale; 227c4762a1bSJed Brown PetscReal ssa_friction_scale; 228c4762a1bSJed Brown PetscReal inertia; 229c4762a1bSJed Brown PRange eta; 230c4762a1bSJed Brown PRange beta2; 231c4762a1bSJed Brown struct { 232c4762a1bSJed Brown PetscReal Bd2,eps,exponent,glen_n; 233c4762a1bSJed Brown } viscosity; 234c4762a1bSJed Brown struct { 235c4762a1bSJed Brown PetscReal irefgam,eps2,exponent; 236c4762a1bSJed Brown } friction; 237c4762a1bSJed Brown struct { 238c4762a1bSJed Brown PetscReal rate,exponent,refvel; 239c4762a1bSJed Brown } erosion; 240c4762a1bSJed Brown PetscReal rhog; 241c4762a1bSJed Brown PetscBool no_slip; 242c4762a1bSJed Brown PetscBool verbose; 243c4762a1bSJed Brown char *mattype; 244c4762a1bSJed Brown char *monitor_basename; 245c4762a1bSJed Brown PetscInt monitor_interval; 246c4762a1bSJed Brown }; 247c4762a1bSJed Brown 248c4762a1bSJed Brown struct _n_Units { 249c4762a1bSJed Brown /* fundamental */ 250c4762a1bSJed Brown PetscReal meter; 251c4762a1bSJed Brown PetscReal kilogram; 252c4762a1bSJed Brown PetscReal second; 253c4762a1bSJed Brown /* derived */ 254c4762a1bSJed Brown PetscReal Pascal; 255c4762a1bSJed Brown PetscReal year; 256c4762a1bSJed Brown }; 257c4762a1bSJed Brown 258c4762a1bSJed Brown static void PrmHexGetZ(const PrmNode pn[],PetscInt k,PetscInt zm,PetscReal zn[]) 259c4762a1bSJed Brown { 260c4762a1bSJed Brown const PetscScalar zm1 = zm-1, 261c4762a1bSJed Brown znl[8] = {pn[0].b + pn[0].h*(PetscScalar)k/zm1, 262c4762a1bSJed Brown pn[1].b + pn[1].h*(PetscScalar)k/zm1, 263c4762a1bSJed Brown pn[2].b + pn[2].h*(PetscScalar)k/zm1, 264c4762a1bSJed Brown pn[3].b + pn[3].h*(PetscScalar)k/zm1, 265c4762a1bSJed Brown pn[0].b + pn[0].h*(PetscScalar)(k+1)/zm1, 266c4762a1bSJed Brown pn[1].b + pn[1].h*(PetscScalar)(k+1)/zm1, 267c4762a1bSJed Brown pn[2].b + pn[2].h*(PetscScalar)(k+1)/zm1, 268c4762a1bSJed Brown pn[3].b + pn[3].h*(PetscScalar)(k+1)/zm1}; 269c4762a1bSJed Brown PetscInt i; 270c4762a1bSJed Brown for (i=0; i<8; i++) zn[i] = PetscRealPart(znl[i]); 271c4762a1bSJed Brown } 272c4762a1bSJed Brown 273c4762a1bSJed Brown /* Compute a gradient of all the 2D fields at four quadrature points. Output for [quadrature_point][direction].field_name */ 274c4762a1bSJed Brown static PetscErrorCode QuadComputeGrad4(const PetscReal dphi[][4][2],PetscReal hx,PetscReal hy,const PrmNode pn[4],PrmNode dp[4][2]) 275c4762a1bSJed Brown { 276c4762a1bSJed Brown PetscInt q,i,f; 277c4762a1bSJed Brown const PetscScalar (*restrict pg)[PRMNODE_SIZE] = (const PetscScalar(*)[PRMNODE_SIZE])pn; /* Get generic array pointers to the node */ 278c4762a1bSJed Brown PetscScalar (*restrict dpg)[2][PRMNODE_SIZE] = (PetscScalar(*)[2][PRMNODE_SIZE])dp; 279c4762a1bSJed Brown 280c4762a1bSJed Brown PetscFunctionBeginUser; 2819566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(dpg,4)); 282c4762a1bSJed Brown for (q=0; q<4; q++) { 283c4762a1bSJed Brown for (i=0; i<4; i++) { 284c4762a1bSJed Brown for (f=0; f<PRMNODE_SIZE; f++) { 285c4762a1bSJed Brown dpg[q][0][f] += dphi[q][i][0]/hx * pg[i][f]; 286c4762a1bSJed Brown dpg[q][1][f] += dphi[q][i][1]/hy * pg[i][f]; 287c4762a1bSJed Brown } 288c4762a1bSJed Brown } 289c4762a1bSJed Brown } 290c4762a1bSJed Brown PetscFunctionReturn(0); 291c4762a1bSJed Brown } 292c4762a1bSJed Brown 293c4762a1bSJed Brown static inline PetscReal StaggeredMidpoint2D(PetscScalar a,PetscScalar b,PetscScalar c,PetscScalar d) 294c4762a1bSJed Brown {return 0.5*PetscRealPart(0.75*a + 0.75*b + 0.25*c + 0.25*d);} 295c4762a1bSJed Brown static inline PetscReal UpwindFlux1D(PetscReal u,PetscReal hL,PetscReal hR) 296c4762a1bSJed Brown {return (u > 0) ? hL*u : hR*u;} 297c4762a1bSJed Brown 298c4762a1bSJed 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), \ 299c4762a1bSJed 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)) 300c4762a1bSJed 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), \ 301c4762a1bSJed 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)) 302c4762a1bSJed 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), \ 303c4762a1bSJed 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)) 304c4762a1bSJed 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), \ 305c4762a1bSJed 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)) 306c4762a1bSJed Brown 307c4762a1bSJed Brown static void PrmNodeGetFaceMeasure(const PrmNode **p,PetscInt i,PetscInt j,PetscScalar h[]) 308c4762a1bSJed Brown { 309c4762a1bSJed Brown /* West */ 310c4762a1bSJed Brown h[0] = StaggeredMidpoint2D(p[i][j].h,p[i-1][j].h,p[i-1][j-1].h,p[i][j-1].h); 311c4762a1bSJed Brown h[1] = StaggeredMidpoint2D(p[i][j].h,p[i-1][j].h,p[i-1][j+1].h,p[i][j+1].h); 312c4762a1bSJed Brown /* East */ 313c4762a1bSJed Brown h[2] = StaggeredMidpoint2D(p[i][j].h,p[i+1][j].h,p[i+1][j+1].h,p[i][j+1].h); 314c4762a1bSJed Brown h[3] = StaggeredMidpoint2D(p[i][j].h,p[i+1][j].h,p[i+1][j-1].h,p[i][j-1].h); 315c4762a1bSJed Brown /* South */ 316c4762a1bSJed Brown h[4] = StaggeredMidpoint2D(p[i][j].h,p[i][j-1].h,p[i+1][j-1].h,p[i+1][j].h); 317c4762a1bSJed Brown h[5] = StaggeredMidpoint2D(p[i][j].h,p[i][j-1].h,p[i-1][j-1].h,p[i-1][j].h); 318c4762a1bSJed Brown /* North */ 319c4762a1bSJed Brown h[6] = StaggeredMidpoint2D(p[i][j].h,p[i][j+1].h,p[i-1][j+1].h,p[i-1][j].h); 320c4762a1bSJed Brown h[7] = StaggeredMidpoint2D(p[i][j].h,p[i][j+1].h,p[i+1][j+1].h,p[i+1][j].h); 321c4762a1bSJed Brown } 322c4762a1bSJed Brown 323c4762a1bSJed Brown /* Tests A and C are from the ISMIP-HOM paper (Pattyn et al. 2008) */ 324c4762a1bSJed Brown static void THIInitialize_HOM_A(THI thi,PetscReal x,PetscReal y,PrmNode *p) 325c4762a1bSJed Brown { 326c4762a1bSJed Brown Units units = thi->units; 327c4762a1bSJed Brown PetscReal s = -x*PetscSinReal(thi->alpha); 328c4762a1bSJed Brown p->b = s - 1000*units->meter + 500*units->meter * PetscSinReal(x*2*PETSC_PI/thi->Lx) * PetscSinReal(y*2*PETSC_PI/thi->Ly); 329c4762a1bSJed Brown p->h = s - p->b; 330c4762a1bSJed 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 */ 331c4762a1bSJed Brown } 332c4762a1bSJed Brown 333c4762a1bSJed Brown static void THIInitialize_HOM_C(THI thi,PetscReal x,PetscReal y,PrmNode *p) 334c4762a1bSJed Brown { 335c4762a1bSJed Brown Units units = thi->units; 336c4762a1bSJed Brown PetscReal s = -x*PetscSinReal(thi->alpha); 337c4762a1bSJed Brown p->b = s - 1000*units->meter; 338c4762a1bSJed Brown p->h = s - p->b; 339c4762a1bSJed Brown /* tau_b = beta2 v is a stress (Pa). 340c4762a1bSJed 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. */ 341c4762a1bSJed 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; 342c4762a1bSJed Brown } 343c4762a1bSJed Brown 344c4762a1bSJed Brown /* These are just toys */ 345c4762a1bSJed Brown 346c4762a1bSJed Brown /* From Fred Herman */ 347c4762a1bSJed Brown static void THIInitialize_HOM_F(THI thi,PetscReal x,PetscReal y,PrmNode *p) 348c4762a1bSJed Brown { 349c4762a1bSJed Brown Units units = thi->units; 350c4762a1bSJed Brown PetscReal s = -x*PetscSinReal(thi->alpha); 351c4762a1bSJed Brown p->b = s - 1000*units->meter + 100*units->meter * PetscSinReal(x*2*PETSC_PI/thi->Lx);/* * sin(y*2*PETSC_PI/thi->Ly); */ 352c4762a1bSJed Brown p->h = s - p->b; 353c4762a1bSJed Brown p->h = (1-(atan((x-thi->Lx/2)/1.)+PETSC_PI/2.)/PETSC_PI)*500*units->meter+1*units->meter; 354c4762a1bSJed Brown s = PetscRealPart(p->b + p->h); 355c4762a1bSJed Brown p->beta2 = -1e-10; 356c4762a1bSJed Brown /* p->beta2 = 1000 * units->Pascal * units->year / units->meter; */ 357c4762a1bSJed Brown } 358c4762a1bSJed Brown 359c4762a1bSJed Brown /* Same bed as test A, free slip everywhere except for a discontinuous jump to a circular sticky region in the middle. */ 360c4762a1bSJed Brown static void THIInitialize_HOM_X(THI thi,PetscReal xx,PetscReal yy,PrmNode *p) 361c4762a1bSJed Brown { 362c4762a1bSJed Brown Units units = thi->units; 363c4762a1bSJed Brown PetscReal x = xx*2*PETSC_PI/thi->Lx - PETSC_PI,y = yy*2*PETSC_PI/thi->Ly - PETSC_PI; /* [-pi,pi] */ 364c4762a1bSJed Brown PetscReal r = PetscSqrtReal(x*x + y*y),s = -x*PetscSinReal(thi->alpha); 365c4762a1bSJed Brown p->b = s - 1000*units->meter + 500*units->meter * PetscSinReal(x + PETSC_PI) * PetscSinReal(y + PETSC_PI); 366c4762a1bSJed Brown p->h = s - p->b; 367c4762a1bSJed Brown p->beta2 = 1000 * (r < 1 ? 2 : 0) * units->Pascal * units->year / units->meter / thi->rhog; 368c4762a1bSJed Brown } 369c4762a1bSJed Brown 370c4762a1bSJed Brown /* Like Z, but with 200 meter cliffs */ 371c4762a1bSJed Brown static void THIInitialize_HOM_Y(THI thi,PetscReal xx,PetscReal yy,PrmNode *p) 372c4762a1bSJed Brown { 373c4762a1bSJed Brown Units units = thi->units; 374c4762a1bSJed Brown PetscReal x = xx*2*PETSC_PI/thi->Lx - PETSC_PI,y = yy*2*PETSC_PI/thi->Ly - PETSC_PI; /* [-pi,pi] */ 375c4762a1bSJed Brown PetscReal r = PetscSqrtReal(x*x + y*y),s = -x*PetscSinReal(thi->alpha); 376c4762a1bSJed Brown p->b = s - 1000*units->meter + 500*units->meter * PetscSinReal(x + PETSC_PI) * PetscSinReal(y + PETSC_PI); 377c4762a1bSJed Brown if (PetscRealPart(p->b) > -700*units->meter) p->b += 200*units->meter; 378c4762a1bSJed Brown p->h = s - p->b; 379c4762a1bSJed 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; 380c4762a1bSJed Brown } 381c4762a1bSJed Brown 382c4762a1bSJed Brown /* Same bed as A, smoothly varying slipperiness, similar to MATLAB's "sombrero" (uncorrelated with bathymetry) */ 383c4762a1bSJed Brown static void THIInitialize_HOM_Z(THI thi,PetscReal xx,PetscReal yy,PrmNode *p) 384c4762a1bSJed Brown { 385c4762a1bSJed Brown Units units = thi->units; 386c4762a1bSJed Brown PetscReal x = xx*2*PETSC_PI/thi->Lx - PETSC_PI,y = yy*2*PETSC_PI/thi->Ly - PETSC_PI; /* [-pi,pi] */ 387c4762a1bSJed Brown PetscReal r = PetscSqrtReal(x*x + y*y),s = -x*PetscSinReal(thi->alpha); 388c4762a1bSJed Brown p->b = s - 1000*units->meter + 500*units->meter * PetscSinReal(x + PETSC_PI) * PetscSinReal(y + PETSC_PI); 389c4762a1bSJed Brown p->h = s - p->b; 390c4762a1bSJed 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; 391c4762a1bSJed Brown } 392c4762a1bSJed Brown 393c4762a1bSJed Brown static void THIFriction(THI thi,PetscReal rbeta2,PetscReal gam,PetscReal *beta2,PetscReal *dbeta2) 394c4762a1bSJed Brown { 395c4762a1bSJed Brown if (thi->friction.irefgam == 0) { 396c4762a1bSJed Brown Units units = thi->units; 397c4762a1bSJed Brown thi->friction.irefgam = 1./(0.5*PetscSqr(100 * units->meter / units->year)); 398c4762a1bSJed Brown thi->friction.eps2 = 0.5*PetscSqr(1.e-4 / thi->friction.irefgam); 399c4762a1bSJed Brown } 400c4762a1bSJed Brown if (thi->friction.exponent == 0) { 401c4762a1bSJed Brown *beta2 = rbeta2; 402c4762a1bSJed Brown *dbeta2 = 0; 403c4762a1bSJed Brown } else { 404c4762a1bSJed Brown *beta2 = rbeta2 * PetscPowReal(thi->friction.eps2 + gam*thi->friction.irefgam,thi->friction.exponent); 405c4762a1bSJed Brown *dbeta2 = thi->friction.exponent * *beta2 / (thi->friction.eps2 + gam*thi->friction.irefgam) * thi->friction.irefgam; 406c4762a1bSJed Brown } 407c4762a1bSJed Brown } 408c4762a1bSJed Brown 409c4762a1bSJed Brown static void THIViscosity(THI thi,PetscReal gam,PetscReal *eta,PetscReal *deta) 410c4762a1bSJed Brown { 411c4762a1bSJed Brown PetscReal Bd2,eps,exponent; 412c4762a1bSJed Brown if (thi->viscosity.Bd2 == 0) { 413c4762a1bSJed Brown Units units = thi->units; 414c4762a1bSJed Brown const PetscReal 415c4762a1bSJed Brown n = thi->viscosity.glen_n, /* Glen exponent */ 416c4762a1bSJed Brown p = 1. + 1./n, /* for Stokes */ 417c4762a1bSJed Brown A = 1.e-16 * PetscPowReal(units->Pascal,-n) / units->year, /* softness parameter (Pa^{-n}/s) */ 418c4762a1bSJed Brown B = PetscPowReal(A,-1./n); /* hardness parameter */ 419c4762a1bSJed Brown thi->viscosity.Bd2 = B/2; 420c4762a1bSJed Brown thi->viscosity.exponent = (p-2)/2; 421c4762a1bSJed Brown thi->viscosity.eps = 0.5*PetscSqr(1e-5 / units->year); 422c4762a1bSJed Brown } 423c4762a1bSJed Brown Bd2 = thi->viscosity.Bd2; 424c4762a1bSJed Brown exponent = thi->viscosity.exponent; 425c4762a1bSJed Brown eps = thi->viscosity.eps; 426c4762a1bSJed Brown *eta = Bd2 * PetscPowReal(eps + gam,exponent); 427c4762a1bSJed Brown *deta = exponent * (*eta) / (eps + gam); 428c4762a1bSJed Brown } 429c4762a1bSJed Brown 430c4762a1bSJed Brown static void THIErosion(THI thi,const Node *vel,PetscScalar *erate,Node *derate) 431c4762a1bSJed Brown { 432c4762a1bSJed Brown const PetscScalar magref2 = 1.e-10 + (PetscSqr(vel->u) + PetscSqr(vel->v)) / PetscSqr(thi->erosion.refvel), 433c4762a1bSJed Brown rate = -thi->erosion.rate*PetscPowScalar(magref2, 0.5*thi->erosion.exponent); 434c4762a1bSJed Brown if (erate) *erate = rate; 435c4762a1bSJed Brown if (derate) { 436c4762a1bSJed Brown if (thi->erosion.exponent == 1) { 437c4762a1bSJed Brown derate->u = 0; 438c4762a1bSJed Brown derate->v = 0; 439c4762a1bSJed Brown } else { 440c4762a1bSJed Brown derate->u = 0.5*thi->erosion.exponent * rate / magref2 * 2. * vel->u / PetscSqr(thi->erosion.refvel); 441c4762a1bSJed Brown derate->v = 0.5*thi->erosion.exponent * rate / magref2 * 2. * vel->v / PetscSqr(thi->erosion.refvel); 442c4762a1bSJed Brown } 443c4762a1bSJed Brown } 444c4762a1bSJed Brown } 445c4762a1bSJed Brown 446c4762a1bSJed Brown static void RangeUpdate(PetscReal *min,PetscReal *max,PetscReal x) 447c4762a1bSJed Brown { 448c4762a1bSJed Brown if (x < *min) *min = x; 449c4762a1bSJed Brown if (x > *max) *max = x; 450c4762a1bSJed Brown } 451c4762a1bSJed Brown 452c4762a1bSJed Brown static void PRangeClear(PRange *p) 453c4762a1bSJed Brown { 454c4762a1bSJed Brown p->cmin = p->min = 1e100; 455c4762a1bSJed Brown p->cmax = p->max = -1e100; 456c4762a1bSJed Brown } 457c4762a1bSJed Brown 458c4762a1bSJed Brown static PetscErrorCode PRangeMinMax(PRange *p,PetscReal min,PetscReal max) 459c4762a1bSJed Brown { 460c4762a1bSJed Brown PetscFunctionBeginUser; 461c4762a1bSJed Brown p->cmin = min; 462c4762a1bSJed Brown p->cmax = max; 463c4762a1bSJed Brown if (min < p->min) p->min = min; 464c4762a1bSJed Brown if (max > p->max) p->max = max; 465c4762a1bSJed Brown PetscFunctionReturn(0); 466c4762a1bSJed Brown } 467c4762a1bSJed Brown 468c4762a1bSJed Brown static PetscErrorCode THIDestroy(THI *thi) 469c4762a1bSJed Brown { 470c4762a1bSJed Brown PetscFunctionBeginUser; 471c4762a1bSJed Brown if (--((PetscObject)(*thi))->refct > 0) PetscFunctionReturn(0); 4729566063dSJacob Faibussowitsch PetscCall(PetscFree((*thi)->units)); 4739566063dSJacob Faibussowitsch PetscCall(PetscFree((*thi)->mattype)); 4749566063dSJacob Faibussowitsch PetscCall(PetscFree((*thi)->monitor_basename)); 4759566063dSJacob Faibussowitsch PetscCall(PetscHeaderDestroy(thi)); 476c4762a1bSJed Brown PetscFunctionReturn(0); 477c4762a1bSJed Brown } 478c4762a1bSJed Brown 479c4762a1bSJed Brown static PetscErrorCode THICreate(MPI_Comm comm,THI *inthi) 480c4762a1bSJed Brown { 481c4762a1bSJed Brown static PetscBool registered = PETSC_FALSE; 482c4762a1bSJed Brown THI thi; 483c4762a1bSJed Brown Units units; 484c4762a1bSJed Brown char monitor_basename[PETSC_MAX_PATH_LEN] = "thi-"; 485c4762a1bSJed Brown PetscErrorCode ierr; 486c4762a1bSJed Brown 487c4762a1bSJed Brown PetscFunctionBeginUser; 488c4762a1bSJed Brown *inthi = 0; 489c4762a1bSJed Brown if (!registered) { 4909566063dSJacob Faibussowitsch PetscCall(PetscClassIdRegister("Toy Hydrostatic Ice",&THI_CLASSID)); 491c4762a1bSJed Brown registered = PETSC_TRUE; 492c4762a1bSJed Brown } 4939566063dSJacob Faibussowitsch PetscCall(PetscHeaderCreate(thi,THI_CLASSID,"THI","Toy Hydrostatic Ice","THI",comm,THIDestroy,0)); 494c4762a1bSJed Brown 4959566063dSJacob Faibussowitsch PetscCall(PetscNew(&thi->units)); 496c4762a1bSJed Brown 497c4762a1bSJed Brown units = thi->units; 498c4762a1bSJed Brown units->meter = 1e-2; 499c4762a1bSJed Brown units->second = 1e-7; 500c4762a1bSJed Brown units->kilogram = 1e-12; 501c4762a1bSJed Brown 502d0609cedSBarry Smith PetscOptionsBegin(comm,NULL,"Scaled units options",""); 503c4762a1bSJed Brown { 5049566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-units_meter","1 meter in scaled length units","",units->meter,&units->meter,NULL)); 5059566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-units_second","1 second in scaled time units","",units->second,&units->second,NULL)); 5069566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-units_kilogram","1 kilogram in scaled mass units","",units->kilogram,&units->kilogram,NULL)); 507c4762a1bSJed Brown } 508d0609cedSBarry Smith PetscOptionsEnd(); 509c4762a1bSJed Brown units->Pascal = units->kilogram / (units->meter * PetscSqr(units->second)); 510c4762a1bSJed Brown units->year = 31556926. * units->second, /* seconds per year */ 511c4762a1bSJed Brown 512c4762a1bSJed Brown thi->Lx = 10.e3; 513c4762a1bSJed Brown thi->Ly = 10.e3; 514c4762a1bSJed Brown thi->Lz = 1000; 515c4762a1bSJed Brown thi->nlevels = 1; 516c4762a1bSJed Brown thi->dirichlet_scale = 1; 517c4762a1bSJed Brown thi->verbose = PETSC_FALSE; 518c4762a1bSJed Brown 519c4762a1bSJed Brown thi->viscosity.glen_n = 3.; 520c4762a1bSJed Brown thi->erosion.rate = 1e-3; /* m/a */ 521c4762a1bSJed Brown thi->erosion.exponent = 1.; 522c4762a1bSJed Brown thi->erosion.refvel = 1.; /* m/a */ 523c4762a1bSJed Brown 524d0609cedSBarry Smith PetscOptionsBegin(comm,NULL,"Toy Hydrostatic Ice options",""); 525c4762a1bSJed Brown { 526c4762a1bSJed Brown QuadratureType quad = QUAD_GAUSS; 527c4762a1bSJed Brown char homexp[] = "A"; 528c4762a1bSJed Brown char mtype[256] = MATSBAIJ; 529c4762a1bSJed Brown PetscReal L,m = 1.0; 530c4762a1bSJed Brown PetscBool flg; 531c4762a1bSJed Brown L = thi->Lx; 5329566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-thi_L","Domain size (m)","",L,&L,&flg)); 533c4762a1bSJed Brown if (flg) thi->Lx = thi->Ly = L; 5349566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-thi_Lx","X Domain size (m)","",thi->Lx,&thi->Lx,NULL)); 5359566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-thi_Ly","Y Domain size (m)","",thi->Ly,&thi->Ly,NULL)); 5369566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-thi_Lz","Z Domain size (m)","",thi->Lz,&thi->Lz,NULL)); 5379566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-thi_hom","ISMIP-HOM experiment (A or C)","",homexp,homexp,sizeof(homexp),NULL)); 538c4762a1bSJed Brown switch (homexp[0] = toupper(homexp[0])) { 539c4762a1bSJed Brown case 'A': 540c4762a1bSJed Brown thi->initialize = THIInitialize_HOM_A; 541c4762a1bSJed Brown thi->no_slip = PETSC_TRUE; 542c4762a1bSJed Brown thi->alpha = 0.5; 543c4762a1bSJed Brown break; 544c4762a1bSJed Brown case 'C': 545c4762a1bSJed Brown thi->initialize = THIInitialize_HOM_C; 546c4762a1bSJed Brown thi->no_slip = PETSC_FALSE; 547c4762a1bSJed Brown thi->alpha = 0.1; 548c4762a1bSJed Brown break; 549c4762a1bSJed Brown case 'F': 550c4762a1bSJed Brown thi->initialize = THIInitialize_HOM_F; 551c4762a1bSJed Brown thi->no_slip = PETSC_FALSE; 552c4762a1bSJed Brown thi->alpha = 0.5; 553c4762a1bSJed Brown break; 554c4762a1bSJed Brown case 'X': 555c4762a1bSJed Brown thi->initialize = THIInitialize_HOM_X; 556c4762a1bSJed Brown thi->no_slip = PETSC_FALSE; 557c4762a1bSJed Brown thi->alpha = 0.3; 558c4762a1bSJed Brown break; 559c4762a1bSJed Brown case 'Y': 560c4762a1bSJed Brown thi->initialize = THIInitialize_HOM_Y; 561c4762a1bSJed Brown thi->no_slip = PETSC_FALSE; 562c4762a1bSJed Brown thi->alpha = 0.5; 563c4762a1bSJed Brown break; 564c4762a1bSJed Brown case 'Z': 565c4762a1bSJed Brown thi->initialize = THIInitialize_HOM_Z; 566c4762a1bSJed Brown thi->no_slip = PETSC_FALSE; 567c4762a1bSJed Brown thi->alpha = 0.5; 568c4762a1bSJed Brown break; 569c4762a1bSJed Brown default: 57098921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"HOM experiment '%c' not implemented",homexp[0]); 571c4762a1bSJed Brown } 5729566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-thi_quadrature","Quadrature to use for 3D elements","",QuadratureTypes,(PetscEnum)quad,(PetscEnum*)&quad,NULL)); 573c4762a1bSJed Brown switch (quad) { 574c4762a1bSJed Brown case QUAD_GAUSS: 575c4762a1bSJed Brown HexQInterp = HexQInterp_Gauss; 576c4762a1bSJed Brown HexQDeriv = HexQDeriv_Gauss; 577c4762a1bSJed Brown break; 578c4762a1bSJed Brown case QUAD_LOBATTO: 579c4762a1bSJed Brown HexQInterp = HexQInterp_Lobatto; 580c4762a1bSJed Brown HexQDeriv = HexQDeriv_Lobatto; 581c4762a1bSJed Brown break; 582c4762a1bSJed Brown } 5839566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-thi_alpha","Bed angle (degrees)","",thi->alpha,&thi->alpha,NULL)); 5849566063dSJacob Faibussowitsch PetscCall(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)); 5859566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-thi_friction_m","Friction exponent, 0=Coulomb, 1=Navier","",m,&m,NULL)); 586c4762a1bSJed Brown thi->friction.exponent = (m-1)/2; 5879566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-thi_erosion_rate","Rate of erosion relative to sliding velocity at reference velocity (m/a)",NULL,thi->erosion.rate,&thi->erosion.rate,NULL)); 5889566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-thi_erosion_exponent","Power of sliding velocity appearing in erosion relation",NULL,thi->erosion.exponent,&thi->erosion.exponent,NULL)); 5899566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-thi_erosion_refvel","Reference sliding velocity for erosion (m/a)",NULL,thi->erosion.refvel,&thi->erosion.refvel,NULL)); 590c4762a1bSJed Brown thi->erosion.rate *= units->meter / units->year; 591c4762a1bSJed Brown thi->erosion.refvel *= units->meter / units->year; 5929566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-thi_dirichlet_scale","Scale Dirichlet boundary conditions by this factor","",thi->dirichlet_scale,&thi->dirichlet_scale,NULL)); 5939566063dSJacob Faibussowitsch PetscCall(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)); 5949566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-thi_inertia","Coefficient of accelaration term in velocity system, physical is almost zero",NULL,thi->inertia,&thi->inertia,NULL)); 5959566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-thi_nlevels","Number of levels of refinement","",thi->nlevels,&thi->nlevels,NULL)); 5969566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-thi_mat_type","Matrix type","MatSetType",MatList,mtype,(char*)mtype,sizeof(mtype),NULL)); 5979566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(mtype,&thi->mattype)); 5989566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-thi_verbose","Enable verbose output (like matrix sizes and statistics)","",thi->verbose,&thi->verbose,NULL)); 5999566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-thi_monitor","Basename to write state files to",NULL,monitor_basename,monitor_basename,sizeof(monitor_basename),&flg)); 600c4762a1bSJed Brown if (flg) { 6019566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(monitor_basename,&thi->monitor_basename)); 602c4762a1bSJed Brown thi->monitor_interval = 1; 6039566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-thi_monitor_interval","Frequency at which to write state files",NULL,thi->monitor_interval,&thi->monitor_interval,NULL)); 604c4762a1bSJed Brown } 605c4762a1bSJed Brown } 606d0609cedSBarry Smith PetscOptionsEnd(); 607c4762a1bSJed Brown 608c4762a1bSJed Brown /* dimensionalize */ 609c4762a1bSJed Brown thi->Lx *= units->meter; 610c4762a1bSJed Brown thi->Ly *= units->meter; 611c4762a1bSJed Brown thi->Lz *= units->meter; 612c4762a1bSJed Brown thi->alpha *= PETSC_PI / 180; 613c4762a1bSJed Brown 614c4762a1bSJed Brown PRangeClear(&thi->eta); 615c4762a1bSJed Brown PRangeClear(&thi->beta2); 616c4762a1bSJed Brown 617c4762a1bSJed Brown { 618c4762a1bSJed Brown PetscReal u = 1000*units->meter/(3e7*units->second), 619c4762a1bSJed Brown gradu = u / (100*units->meter),eta,deta, 620c4762a1bSJed Brown rho = 910 * units->kilogram/PetscPowRealInt(units->meter,3), 621c4762a1bSJed Brown grav = 9.81 * units->meter/PetscSqr(units->second), 622c4762a1bSJed Brown driving = rho * grav * PetscSinReal(thi->alpha) * 1000*units->meter; 623c4762a1bSJed Brown THIViscosity(thi,0.5*gradu*gradu,&eta,&deta); 624c4762a1bSJed Brown thi->rhog = rho * grav; 625c4762a1bSJed Brown if (thi->verbose) { 62663a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)thi),"Units: meter %8.2g second %8.2g kg %8.2g Pa %8.2g\n",(double)units->meter,(double)units->second,(double)units->kilogram,(double)units->Pascal)); 62763a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)thi),"Domain (%6.2g,%6.2g,%6.2g), pressure %8.2g, driving stress %8.2g\n",(double)thi->Lx,(double)thi->Ly,(double)thi->Lz,(double)(rho*grav*1e3*units->meter),(double)driving)); 62863a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)thi),"Large velocity 1km/a %8.2g, velocity gradient %8.2g, eta %8.2g, stress %8.2g, ratio %8.2g\n",(double)u,(double)gradu,(double)eta,(double)(2*eta*gradu,2*eta*gradu/driving))); 629c4762a1bSJed Brown THIViscosity(thi,0.5*PetscSqr(1e-3*gradu),&eta,&deta); 63063a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)thi),"Small velocity 1m/a %8.2g, velocity gradient %8.2g, eta %8.2g, stress %8.2g, ratio %8.2g\n",(double)(1e-3*u),(double)(1e-3*gradu),(double)eta,(double)(2*eta*1e-3*gradu,2*eta*1e-3*gradu/driving))); 631c4762a1bSJed Brown } 632c4762a1bSJed Brown } 633c4762a1bSJed Brown 634c4762a1bSJed Brown *inthi = thi; 635c4762a1bSJed Brown PetscFunctionReturn(0); 636c4762a1bSJed Brown } 637c4762a1bSJed Brown 638c4762a1bSJed Brown /* Our problem is periodic, but the domain has a mean slope of alpha so the bed does not line up between the upstream 639c4762a1bSJed Brown * and downstream ends of the domain. This function fixes the ghost values so that the domain appears truly periodic in 640c4762a1bSJed Brown * the horizontal. */ 641c4762a1bSJed Brown static PetscErrorCode THIFixGhosts(THI thi,DM da3,DM da2,Vec X3,Vec X2) 642c4762a1bSJed Brown { 643c4762a1bSJed Brown DMDALocalInfo info; 644c4762a1bSJed Brown PrmNode **x2; 645c4762a1bSJed Brown PetscInt i,j; 646c4762a1bSJed Brown 647c4762a1bSJed Brown PetscFunctionBeginUser; 6489566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da3,&info)); 6499566063dSJacob Faibussowitsch /* PetscCall(VecView(X2,PETSC_VIEWER_STDOUT_WORLD)); */ 6509566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da2,X2,&x2)); 651c4762a1bSJed Brown for (i=info.gzs; i<info.gzs+info.gzm; i++) { 652c4762a1bSJed Brown if (i > -1 && i < info.mz) continue; 653c4762a1bSJed Brown for (j=info.gys; j<info.gys+info.gym; j++) { 654c4762a1bSJed Brown x2[i][j].b += PetscSinReal(thi->alpha) * thi->Lx * (i<0 ? 1.0 : -1.0); 655c4762a1bSJed Brown } 656c4762a1bSJed Brown } 6579566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da2,X2,&x2)); 6589566063dSJacob Faibussowitsch /* PetscCall(VecView(X2,PETSC_VIEWER_STDOUT_WORLD)); */ 659c4762a1bSJed Brown PetscFunctionReturn(0); 660c4762a1bSJed Brown } 661c4762a1bSJed Brown 662c4762a1bSJed Brown static PetscErrorCode THIInitializePrm(THI thi,DM da2prm,PrmNode **p) 663c4762a1bSJed Brown { 664c4762a1bSJed Brown PetscInt i,j,xs,xm,ys,ym,mx,my; 665c4762a1bSJed Brown 666c4762a1bSJed Brown PetscFunctionBeginUser; 6679566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da2prm,&ys,&xs,0,&ym,&xm,0)); 6689566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da2prm,0, &my,&mx,0, 0,0,0, 0,0,0,0,0,0)); 669c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 670c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 671c4762a1bSJed Brown PetscReal xx = thi->Lx*i/mx,yy = thi->Ly*j/my; 672c4762a1bSJed Brown thi->initialize(thi,xx,yy,&p[i][j]); 673c4762a1bSJed Brown } 674c4762a1bSJed Brown } 675c4762a1bSJed Brown PetscFunctionReturn(0); 676c4762a1bSJed Brown } 677c4762a1bSJed Brown 678c4762a1bSJed Brown static PetscErrorCode THIInitial(THI thi,DM pack,Vec X) 679c4762a1bSJed Brown { 680c4762a1bSJed Brown DM da3,da2; 681c4762a1bSJed Brown PetscInt i,j,k,xs,xm,ys,ym,zs,zm,mx,my; 682c4762a1bSJed Brown PetscReal hx,hy; 683c4762a1bSJed Brown PrmNode **prm; 684c4762a1bSJed Brown Node ***x; 685c4762a1bSJed Brown Vec X3g,X2g,X2; 686c4762a1bSJed Brown 687c4762a1bSJed Brown PetscFunctionBeginUser; 6889566063dSJacob Faibussowitsch PetscCall(DMCompositeGetEntries(pack,&da3,&da2)); 6899566063dSJacob Faibussowitsch PetscCall(DMCompositeGetAccess(pack,X,&X3g,&X2g)); 6909566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da2,&X2)); 691c4762a1bSJed Brown 6929566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da3,0, 0,&my,&mx, 0,0,0, 0,0,0,0,0,0)); 6939566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da3,&zs,&ys,&xs,&zm,&ym,&xm)); 6949566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da3,X3g,&x)); 6959566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da2,X2,&prm)); 696c4762a1bSJed Brown 6979566063dSJacob Faibussowitsch PetscCall(THIInitializePrm(thi,da2,prm)); 698c4762a1bSJed Brown 699c4762a1bSJed Brown hx = thi->Lx / mx; 700c4762a1bSJed Brown hy = thi->Ly / my; 701c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 702c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 703c4762a1bSJed Brown for (k=zs; k<zs+zm; k++) { 704c4762a1bSJed Brown const PetscScalar zm1 = zm-1, 705c4762a1bSJed 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), 706c4762a1bSJed 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); 707c4762a1bSJed Brown x[i][j][k].u = 0. * drivingx * prm[i][j].h*(PetscScalar)k/zm1; 708c4762a1bSJed Brown x[i][j][k].v = 0. * drivingy * prm[i][j].h*(PetscScalar)k/zm1; 709c4762a1bSJed Brown } 710c4762a1bSJed Brown } 711c4762a1bSJed Brown } 712c4762a1bSJed Brown 7139566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da3,X3g,&x)); 7149566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da2,X2,&prm)); 715c4762a1bSJed Brown 7169566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(da2,X2,INSERT_VALUES,X2g)); 7179566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd (da2,X2,INSERT_VALUES,X2g)); 7189566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da2,&X2)); 719c4762a1bSJed Brown 7209566063dSJacob Faibussowitsch PetscCall(DMCompositeRestoreAccess(pack,X,&X3g,&X2g)); 721c4762a1bSJed Brown PetscFunctionReturn(0); 722c4762a1bSJed Brown } 723c4762a1bSJed Brown 724c4762a1bSJed 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) 725c4762a1bSJed Brown { 726c4762a1bSJed Brown PetscInt l,ll; 727c4762a1bSJed Brown PetscScalar gam; 728c4762a1bSJed Brown 729c4762a1bSJed Brown du[0] = du[1] = du[2] = 0; 730c4762a1bSJed Brown dv[0] = dv[1] = dv[2] = 0; 731c4762a1bSJed Brown *u = 0; 732c4762a1bSJed Brown *v = 0; 733c4762a1bSJed Brown for (l=0; l<8; l++) { 734c4762a1bSJed Brown *u += phi[l] * n[l].u; 735c4762a1bSJed Brown *v += phi[l] * n[l].v; 736c4762a1bSJed Brown for (ll=0; ll<3; ll++) { 737c4762a1bSJed Brown du[ll] += dphi[l][ll] * n[l].u; 738c4762a1bSJed Brown dv[ll] += dphi[l][ll] * n[l].v; 739c4762a1bSJed Brown } 740c4762a1bSJed Brown } 741c4762a1bSJed 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]); 742c4762a1bSJed Brown THIViscosity(thi,PetscRealPart(gam),eta,deta); 743c4762a1bSJed Brown } 744c4762a1bSJed Brown 745c4762a1bSJed Brown static PetscErrorCode THIFunctionLocal_3D(DMDALocalInfo *info,const Node ***x,const PrmNode **prm,const Node ***xdot,Node ***f,THI thi) 746c4762a1bSJed Brown { 747c4762a1bSJed Brown PetscInt xs,ys,xm,ym,zm,i,j,k,q,l; 748c4762a1bSJed Brown PetscReal hx,hy,etamin,etamax,beta2min,beta2max; 749c4762a1bSJed Brown 750c4762a1bSJed Brown PetscFunctionBeginUser; 751c4762a1bSJed Brown xs = info->zs; 752c4762a1bSJed Brown ys = info->ys; 753c4762a1bSJed Brown xm = info->zm; 754c4762a1bSJed Brown ym = info->ym; 755c4762a1bSJed Brown zm = info->xm; 756c4762a1bSJed Brown hx = thi->Lx / info->mz; 757c4762a1bSJed Brown hy = thi->Ly / info->my; 758c4762a1bSJed Brown 759c4762a1bSJed Brown etamin = 1e100; 760c4762a1bSJed Brown etamax = 0; 761c4762a1bSJed Brown beta2min = 1e100; 762c4762a1bSJed Brown beta2max = 0; 763c4762a1bSJed Brown 764c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 765c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 766c4762a1bSJed Brown PrmNode pn[4],dpn[4][2]; 767c4762a1bSJed Brown QuadExtract(prm,i,j,pn); 7689566063dSJacob Faibussowitsch PetscCall(QuadComputeGrad4(QuadQDeriv,hx,hy,pn,dpn)); 769c4762a1bSJed Brown for (k=0; k<zm-1; k++) { 770c4762a1bSJed Brown PetscInt ls = 0; 771c4762a1bSJed Brown Node n[8],ndot[8],*fn[8]; 772c4762a1bSJed Brown PetscReal zn[8],etabase = 0; 7732f613bf5SBarry Smith 774c4762a1bSJed Brown PrmHexGetZ(pn,k,zm,zn); 775c4762a1bSJed Brown HexExtract(x,i,j,k,n); 7762f613bf5SBarry Smith HexExtract(xdot,i,j,k,ndot); 777c4762a1bSJed Brown HexExtractRef(f,i,j,k,fn); 778c4762a1bSJed Brown if (thi->no_slip && k == 0) { 779c4762a1bSJed Brown for (l=0; l<4; l++) n[l].u = n[l].v = 0; 780c4762a1bSJed Brown /* The first 4 basis functions lie on the bottom layer, so their contribution is exactly 0, hence we can skip them */ 781c4762a1bSJed Brown ls = 4; 782c4762a1bSJed Brown } 783c4762a1bSJed Brown for (q=0; q<8; q++) { 784c4762a1bSJed Brown PetscReal dz[3],phi[8],dphi[8][3],jw,eta,deta; 785c4762a1bSJed Brown PetscScalar du[3],dv[3],u,v,udot=0,vdot=0; 786c4762a1bSJed Brown for (l=ls; l<8; l++) { 787c4762a1bSJed Brown udot += HexQInterp[q][l]*ndot[l].u; 788c4762a1bSJed Brown vdot += HexQInterp[q][l]*ndot[l].v; 789c4762a1bSJed Brown } 790c4762a1bSJed Brown HexGrad(HexQDeriv[q],zn,dz); 791c4762a1bSJed Brown HexComputeGeometry(q,hx,hy,dz,phi,dphi,&jw); 792c4762a1bSJed Brown PointwiseNonlinearity(thi,n,phi,dphi,&u,&v,du,dv,&eta,&deta); 793c4762a1bSJed Brown jw /= thi->rhog; /* scales residuals to be O(1) */ 794c4762a1bSJed Brown if (q == 0) etabase = eta; 795c4762a1bSJed Brown RangeUpdate(&etamin,&etamax,eta); 796c4762a1bSJed Brown for (l=ls; l<8; l++) { /* test functions */ 797c4762a1bSJed 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}; 798c4762a1bSJed Brown const PetscReal pp = phi[l],*dp = dphi[l]; 799c4762a1bSJed 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]; 800c4762a1bSJed 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]; 801c4762a1bSJed Brown fn[l]->u += pp*jw*udot*thi->inertia*pp; 802c4762a1bSJed Brown fn[l]->v += pp*jw*vdot*thi->inertia*pp; 803c4762a1bSJed Brown } 804c4762a1bSJed Brown } 805c4762a1bSJed Brown if (k == 0) { /* we are on a bottom face */ 806c4762a1bSJed Brown if (thi->no_slip) { 807c4762a1bSJed Brown /* Note: Non-Galerkin coarse grid operators are very sensitive to the scaling of Dirichlet boundary 808c4762a1bSJed Brown * conditions. After shenanigans above, etabase contains the effective viscosity at the closest quadrature 809c4762a1bSJed Brown * point to the bed. We want the diagonal entry in the Dirichlet condition to have similar magnitude to the 810c4762a1bSJed Brown * diagonal entry corresponding to the adjacent node. The fundamental scaling of the viscous part is in 811c4762a1bSJed Brown * diagu, diagv below. This scaling is easy to recognize by considering the finite difference operator after 812c4762a1bSJed Brown * scaling by element size. The no-slip Dirichlet condition is scaled by this factor, and also in the 813c4762a1bSJed Brown * assembled matrix (see the similar block in THIJacobianLocal). 814c4762a1bSJed Brown * 815c4762a1bSJed Brown * Note that the residual at this Dirichlet node is linear in the state at this node, but also depends 816c4762a1bSJed Brown * (nonlinearly in general) on the neighboring interior nodes through the local viscosity. This will make 817c4762a1bSJed Brown * a matrix-free Jacobian have extra entries in the corresponding row. We assemble only the diagonal part, 818c4762a1bSJed Brown * so the solution will exactly satisfy the boundary condition after the first linear iteration. 819c4762a1bSJed Brown */ 820c4762a1bSJed Brown const PetscReal hz = PetscRealPart(pn[0].h)/(zm-1.); 821c4762a1bSJed 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); 822c4762a1bSJed Brown fn[0]->u = thi->dirichlet_scale*diagu*x[i][j][k].u; 823c4762a1bSJed Brown fn[0]->v = thi->dirichlet_scale*diagv*x[i][j][k].v; 824c4762a1bSJed Brown } else { /* Integrate over bottom face to apply boundary condition */ 825c4762a1bSJed 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) */ 826c4762a1bSJed Brown const PetscReal jw = 0.25*hx*hy,*phi = QuadQInterp[q]; 827c4762a1bSJed Brown PetscScalar u =0,v=0,rbeta2=0; 828c4762a1bSJed Brown PetscReal beta2,dbeta2; 829c4762a1bSJed Brown for (l=0; l<4; l++) { 830c4762a1bSJed Brown u += phi[l]*n[l].u; 831c4762a1bSJed Brown v += phi[l]*n[l].v; 832c4762a1bSJed Brown rbeta2 += phi[l]*pn[l].beta2; 833c4762a1bSJed Brown } 834c4762a1bSJed Brown THIFriction(thi,PetscRealPart(rbeta2),PetscRealPart(u*u+v*v)/2,&beta2,&dbeta2); 835c4762a1bSJed Brown RangeUpdate(&beta2min,&beta2max,beta2); 836c4762a1bSJed Brown for (l=0; l<4; l++) { 837c4762a1bSJed Brown const PetscReal pp = phi[l]; 838c4762a1bSJed Brown fn[ls+l]->u += pp*jw*beta2*u; 839c4762a1bSJed Brown fn[ls+l]->v += pp*jw*beta2*v; 840c4762a1bSJed Brown } 841c4762a1bSJed Brown } 842c4762a1bSJed Brown } 843c4762a1bSJed Brown } 844c4762a1bSJed Brown } 845c4762a1bSJed Brown } 846c4762a1bSJed Brown } 847c4762a1bSJed Brown 8489566063dSJacob Faibussowitsch PetscCall(PRangeMinMax(&thi->eta,etamin,etamax)); 8499566063dSJacob Faibussowitsch PetscCall(PRangeMinMax(&thi->beta2,beta2min,beta2max)); 850c4762a1bSJed Brown PetscFunctionReturn(0); 851c4762a1bSJed Brown } 852c4762a1bSJed Brown 853c4762a1bSJed Brown static PetscErrorCode THIFunctionLocal_2D(DMDALocalInfo *info,const Node ***x,const PrmNode **prm,const PrmNode **prmdot,PrmNode **f,THI thi) 854c4762a1bSJed Brown { 855c4762a1bSJed Brown PetscInt xs,ys,xm,ym,zm,i,j,k; 856c4762a1bSJed Brown 857c4762a1bSJed Brown PetscFunctionBeginUser; 858c4762a1bSJed Brown xs = info->zs; 859c4762a1bSJed Brown ys = info->ys; 860c4762a1bSJed Brown xm = info->zm; 861c4762a1bSJed Brown ym = info->ym; 862c4762a1bSJed Brown zm = info->xm; 863c4762a1bSJed Brown 864c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 865c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 866c4762a1bSJed Brown PetscScalar div = 0,erate,h[8]; 867c4762a1bSJed Brown PrmNodeGetFaceMeasure(prm,i,j,h); 868c4762a1bSJed Brown for (k=0; k<zm; k++) { 869c4762a1bSJed Brown PetscScalar weight = (k==0 || k == zm-1) ? 0.5/(zm-1) : 1.0/(zm-1); 870c4762a1bSJed Brown if (0) { /* centered flux */ 871c4762a1bSJed 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) 872c4762a1bSJed 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) 873c4762a1bSJed 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) 874c4762a1bSJed 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) 875c4762a1bSJed 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) 876c4762a1bSJed 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) 877c4762a1bSJed 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) 878c4762a1bSJed 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)); 879c4762a1bSJed Brown } else { /* Upwind flux */ 880c4762a1bSJed Brown div += weight*(-UpwindFluxXW(x,prm,h,i,j,k, 1) 881c4762a1bSJed Brown -UpwindFluxXW(x,prm,h,i,j,k,-1) 882c4762a1bSJed Brown +UpwindFluxXE(x,prm,h,i,j,k, 1) 883c4762a1bSJed Brown +UpwindFluxXE(x,prm,h,i,j,k,-1) 884c4762a1bSJed Brown -UpwindFluxYS(x,prm,h,i,j,k, 1) 885c4762a1bSJed Brown -UpwindFluxYS(x,prm,h,i,j,k,-1) 886c4762a1bSJed Brown +UpwindFluxYN(x,prm,h,i,j,k, 1) 887c4762a1bSJed Brown +UpwindFluxYN(x,prm,h,i,j,k,-1)); 888c4762a1bSJed Brown } 889c4762a1bSJed Brown } 890c4762a1bSJed Brown /* printf("div[%d][%d] %g\n",i,j,div); */ 891c4762a1bSJed Brown THIErosion(thi,&x[i][j][0],&erate,NULL); 892c4762a1bSJed Brown f[i][j].b = prmdot[i][j].b - erate; 893c4762a1bSJed Brown f[i][j].h = prmdot[i][j].h + div; 894c4762a1bSJed Brown f[i][j].beta2 = prmdot[i][j].beta2; 895c4762a1bSJed Brown } 896c4762a1bSJed Brown } 897c4762a1bSJed Brown PetscFunctionReturn(0); 898c4762a1bSJed Brown } 899c4762a1bSJed Brown 900c4762a1bSJed Brown static PetscErrorCode THIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx) 901c4762a1bSJed Brown { 902c4762a1bSJed Brown THI thi = (THI)ctx; 903c4762a1bSJed Brown DM pack,da3,da2; 904c4762a1bSJed Brown Vec X3,X2,Xdot3,Xdot2,F3,F2,F3g,F2g; 905c4762a1bSJed Brown const Node ***x3,***xdot3; 906c4762a1bSJed Brown const PrmNode **x2,**xdot2; 907c4762a1bSJed Brown Node ***f3; 908c4762a1bSJed Brown PrmNode **f2; 909c4762a1bSJed Brown DMDALocalInfo info3; 910c4762a1bSJed Brown 911c4762a1bSJed Brown PetscFunctionBeginUser; 9129566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts,&pack)); 9139566063dSJacob Faibussowitsch PetscCall(DMCompositeGetEntries(pack,&da3,&da2)); 9149566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da3,&info3)); 9159566063dSJacob Faibussowitsch PetscCall(DMCompositeGetLocalVectors(pack,&X3,&X2)); 9169566063dSJacob Faibussowitsch PetscCall(DMCompositeGetLocalVectors(pack,&Xdot3,&Xdot2)); 9179566063dSJacob Faibussowitsch PetscCall(DMCompositeScatter(pack,X,X3,X2)); 9189566063dSJacob Faibussowitsch PetscCall(THIFixGhosts(thi,da3,da2,X3,X2)); 9199566063dSJacob Faibussowitsch PetscCall(DMCompositeScatter(pack,Xdot,Xdot3,Xdot2)); 920c4762a1bSJed Brown 9219566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da3,&F3)); 9229566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da2,&F2)); 9239566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F3)); 924c4762a1bSJed Brown 9259566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da3,X3,&x3)); 9269566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da2,X2,&x2)); 9279566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da3,Xdot3,&xdot3)); 9289566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da2,Xdot2,&xdot2)); 9299566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da3,F3,&f3)); 9309566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da2,F2,&f2)); 931c4762a1bSJed Brown 9329566063dSJacob Faibussowitsch PetscCall(THIFunctionLocal_3D(&info3,x3,x2,xdot3,f3,thi)); 9339566063dSJacob Faibussowitsch PetscCall(THIFunctionLocal_2D(&info3,x3,x2,xdot2,f2,thi)); 934c4762a1bSJed Brown 9359566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da3,X3,&x3)); 9369566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da2,X2,&x2)); 9379566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da3,Xdot3,&xdot3)); 9389566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da2,Xdot2,&xdot2)); 9399566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da3,F3,&f3)); 9409566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da2,F2,&f2)); 941c4762a1bSJed Brown 9429566063dSJacob Faibussowitsch PetscCall(DMCompositeRestoreLocalVectors(pack,&X3,&X2)); 9439566063dSJacob Faibussowitsch PetscCall(DMCompositeRestoreLocalVectors(pack,&Xdot3,&Xdot2)); 944c4762a1bSJed Brown 9459566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F)); 9469566063dSJacob Faibussowitsch PetscCall(DMCompositeGetAccess(pack,F,&F3g,&F2g)); 9479566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(da3,F3,ADD_VALUES,F3g)); 9489566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd (da3,F3,ADD_VALUES,F3g)); 9499566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(da2,F2,INSERT_VALUES,F2g)); 9509566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd (da2,F2,INSERT_VALUES,F2g)); 951c4762a1bSJed Brown 952c4762a1bSJed Brown if (thi->verbose) { 953c4762a1bSJed Brown PetscViewer viewer; 9549566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)thi),&viewer)); 9559566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"3D_Velocity residual (bs=2):\n")); 9569566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 9579566063dSJacob Faibussowitsch PetscCall(VecView(F3,viewer)); 9589566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 9599566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"2D_Fields residual (bs=3):\n")); 9609566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 9619566063dSJacob Faibussowitsch PetscCall(VecView(F2,viewer)); 9629566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 963c4762a1bSJed Brown } 964c4762a1bSJed Brown 9659566063dSJacob Faibussowitsch PetscCall(DMCompositeRestoreAccess(pack,F,&F3g,&F2g)); 966c4762a1bSJed Brown 9679566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da3,&F3)); 9689566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da2,&F2)); 969c4762a1bSJed Brown PetscFunctionReturn(0); 970c4762a1bSJed Brown } 971c4762a1bSJed Brown 972c4762a1bSJed Brown static PetscErrorCode THIMatrixStatistics(THI thi,Mat B,PetscViewer viewer) 973c4762a1bSJed Brown { 974c4762a1bSJed Brown PetscReal nrm; 975c4762a1bSJed Brown PetscInt m; 976c4762a1bSJed Brown PetscMPIInt rank; 977c4762a1bSJed Brown 978c4762a1bSJed Brown PetscFunctionBeginUser; 9799566063dSJacob Faibussowitsch PetscCall(MatNorm(B,NORM_FROBENIUS,&nrm)); 9809566063dSJacob Faibussowitsch PetscCall(MatGetSize(B,&m,0)); 9819566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)B),&rank)); 982dd400576SPatrick Sanan if (rank == 0) { 983c4762a1bSJed Brown PetscScalar val0,val2; 9849566063dSJacob Faibussowitsch PetscCall(MatGetValue(B,0,0,&val0)); 9859566063dSJacob Faibussowitsch PetscCall(MatGetValue(B,2,2,&val2)); 98663a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Matrix dim %8" PetscInt_FMT " norm %8.2e, (0,0) %8.2e (2,2) %8.2e, eta [%8.2e,%8.2e] beta2 [%8.2e,%8.2e]\n",m,(double)nrm,(double)PetscRealPart(val0),(double)PetscRealPart(val2),(double)thi->eta.cmin,(double)thi->eta.cmax,(double)thi->beta2.cmin,(double)thi->beta2.cmax)); 987c4762a1bSJed Brown } 988c4762a1bSJed Brown PetscFunctionReturn(0); 989c4762a1bSJed Brown } 990c4762a1bSJed Brown 991c4762a1bSJed Brown static PetscErrorCode THISurfaceStatistics(DM pack,Vec X,PetscReal *min,PetscReal *max,PetscReal *mean) 992c4762a1bSJed Brown { 993c4762a1bSJed Brown DM da3,da2; 994c4762a1bSJed Brown Vec X3,X2; 995c4762a1bSJed Brown Node ***x; 996c4762a1bSJed Brown PetscInt i,j,xs,ys,zs,xm,ym,zm,mx,my,mz; 997c4762a1bSJed Brown PetscReal umin = 1e100,umax=-1e100; 998c4762a1bSJed Brown PetscScalar usum =0.0,gusum; 999c4762a1bSJed Brown 1000c4762a1bSJed Brown PetscFunctionBeginUser; 10019566063dSJacob Faibussowitsch PetscCall(DMCompositeGetEntries(pack,&da3,&da2)); 10029566063dSJacob Faibussowitsch PetscCall(DMCompositeGetAccess(pack,X,&X3,&X2)); 1003c4762a1bSJed Brown *min = *max = *mean = 0; 10049566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da3,0, &mz,&my,&mx, 0,0,0, 0,0,0,0,0,0)); 10059566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da3,&zs,&ys,&xs,&zm,&ym,&xm)); 10063c633725SBarry Smith PetscCheck(zs == 0 && zm == mz,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Unexpected decomposition"); 10079566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da3,X3,&x)); 1008c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 1009c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 1010c4762a1bSJed Brown PetscReal u = PetscRealPart(x[i][j][zm-1].u); 1011c4762a1bSJed Brown RangeUpdate(&umin,&umax,u); 1012c4762a1bSJed Brown usum += u; 1013c4762a1bSJed Brown } 1014c4762a1bSJed Brown } 10159566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da3,X3,&x)); 10169566063dSJacob Faibussowitsch PetscCall(DMCompositeRestoreAccess(pack,X,&X3,&X2)); 1017c4762a1bSJed Brown 10189566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(&umin,min,1,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)da3))); 10199566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(&umax,max,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)da3))); 10209566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(&usum,&gusum,1,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)da3))); 1021c4762a1bSJed Brown *mean = PetscRealPart(gusum) / (mx*my); 1022c4762a1bSJed Brown PetscFunctionReturn(0); 1023c4762a1bSJed Brown } 1024c4762a1bSJed Brown 1025c4762a1bSJed Brown static PetscErrorCode THISolveStatistics(THI thi,TS ts,PetscInt coarsened,const char name[]) 1026c4762a1bSJed Brown { 1027c4762a1bSJed Brown MPI_Comm comm; 1028c4762a1bSJed Brown DM pack; 1029c4762a1bSJed Brown Vec X,X3,X2; 1030c4762a1bSJed Brown 1031c4762a1bSJed Brown PetscFunctionBeginUser; 10329566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)thi,&comm)); 10339566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts,&pack)); 10349566063dSJacob Faibussowitsch PetscCall(TSGetSolution(ts,&X)); 10359566063dSJacob Faibussowitsch PetscCall(DMCompositeGetAccess(pack,X,&X3,&X2)); 10369566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm,"Solution statistics after solve: %s\n",name)); 1037c4762a1bSJed Brown { 1038c4762a1bSJed Brown PetscInt its,lits; 1039c4762a1bSJed Brown SNESConvergedReason reason; 1040c4762a1bSJed Brown SNES snes; 10419566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts,&snes)); 10429566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes,&its)); 10439566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes,&reason)); 10449566063dSJacob Faibussowitsch PetscCall(SNESGetLinearSolveIterations(snes,&lits)); 104563a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm,"%s: Number of SNES iterations = %" PetscInt_FMT ", total linear iterations = %" PetscInt_FMT "\n",SNESConvergedReasons[reason],its,lits)); 1046c4762a1bSJed Brown } 1047c4762a1bSJed Brown { 1048c4762a1bSJed Brown PetscReal nrm2,tmin[3]={1e100,1e100,1e100},tmax[3]={-1e100,-1e100,-1e100},min[3],max[3]; 1049c4762a1bSJed Brown PetscInt i,j,m; 1050c4762a1bSJed Brown PetscScalar *x; 10519566063dSJacob Faibussowitsch PetscCall(VecNorm(X3,NORM_2,&nrm2)); 10529566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(X3,&m)); 10539566063dSJacob Faibussowitsch PetscCall(VecGetArray(X3,&x)); 1054c4762a1bSJed Brown for (i=0; i<m; i+=2) { 1055c4762a1bSJed Brown PetscReal u = PetscRealPart(x[i]),v = PetscRealPart(x[i+1]),c = PetscSqrtReal(u*u+v*v); 1056c4762a1bSJed Brown tmin[0] = PetscMin(u,tmin[0]); 1057c4762a1bSJed Brown tmin[1] = PetscMin(v,tmin[1]); 1058c4762a1bSJed Brown tmin[2] = PetscMin(c,tmin[2]); 1059c4762a1bSJed Brown tmax[0] = PetscMax(u,tmax[0]); 1060c4762a1bSJed Brown tmax[1] = PetscMax(v,tmax[1]); 1061c4762a1bSJed Brown tmax[2] = PetscMax(c,tmax[2]); 1062c4762a1bSJed Brown } 10639566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(X,&x)); 10649566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(tmin,min,3,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)thi))); 10659566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(tmax,max,3,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)thi))); 1066c4762a1bSJed Brown /* Dimensionalize to meters/year */ 1067c4762a1bSJed Brown nrm2 *= thi->units->year / thi->units->meter; 1068c4762a1bSJed Brown for (j=0; j<3; j++) { 1069c4762a1bSJed Brown min[j] *= thi->units->year / thi->units->meter; 1070c4762a1bSJed Brown max[j] *= thi->units->year / thi->units->meter; 1071c4762a1bSJed Brown } 107263a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm,"|X|_2 %g u in [%g, %g] v in [%g, %g] c in [%g, %g] \n",(double)nrm2,(double)min[0],(double)max[0],(double)min[1],(double)max[1],(double)min[2],(double)max[2])); 1073c4762a1bSJed Brown { 1074c4762a1bSJed Brown PetscReal umin,umax,umean; 10759566063dSJacob Faibussowitsch PetscCall(THISurfaceStatistics(pack,X,&umin,&umax,&umean)); 1076c4762a1bSJed Brown umin *= thi->units->year / thi->units->meter; 1077c4762a1bSJed Brown umax *= thi->units->year / thi->units->meter; 1078c4762a1bSJed Brown umean *= thi->units->year / thi->units->meter; 107963a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm,"Surface statistics: u in [%12.6e, %12.6e] mean %12.6e\n",(double)umin,(double)umax,(double)umean)); 1080c4762a1bSJed Brown } 1081c4762a1bSJed Brown /* These values stay nondimensional */ 108263a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm,"Global eta range [%g, %g], converged range [%g, %g]\n",(double)thi->eta.min,(double)thi->eta.max,(double)thi->eta.cmin,(double)thi->eta.cmax)); 108363a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm,"Global beta2 range [%g, %g], converged range [%g, %g]\n",(double)thi->beta2.min,(double)thi->beta2.max,(double)thi->beta2.cmin,(double)thi->beta2.cmax)); 1084c4762a1bSJed Brown } 10859566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm,"\n")); 10869566063dSJacob Faibussowitsch PetscCall(DMCompositeRestoreAccess(pack,X,&X3,&X2)); 1087c4762a1bSJed Brown PetscFunctionReturn(0); 1088c4762a1bSJed Brown } 1089c4762a1bSJed Brown 1090c4762a1bSJed Brown static inline PetscInt DMDALocalIndex3D(DMDALocalInfo *info,PetscInt i,PetscInt j,PetscInt k) 1091c4762a1bSJed Brown {return ((i-info->gzs)*info->gym + (j-info->gys))*info->gxm + (k-info->gxs);} 1092c4762a1bSJed Brown static inline PetscInt DMDALocalIndex2D(DMDALocalInfo *info,PetscInt i,PetscInt j) 1093c4762a1bSJed Brown {return (i-info->gzs)*info->gym + (j-info->gys);} 1094c4762a1bSJed Brown 1095c4762a1bSJed Brown static PetscErrorCode THIJacobianLocal_Momentum(DMDALocalInfo *info,const Node ***x,const PrmNode **prm,Mat B,Mat Bcpl,THI thi) 1096c4762a1bSJed Brown { 1097c4762a1bSJed Brown PetscInt xs,ys,xm,ym,zm,i,j,k,q,l,ll; 1098c4762a1bSJed Brown PetscReal hx,hy; 1099c4762a1bSJed Brown 1100c4762a1bSJed Brown PetscFunctionBeginUser; 1101c4762a1bSJed Brown xs = info->zs; 1102c4762a1bSJed Brown ys = info->ys; 1103c4762a1bSJed Brown xm = info->zm; 1104c4762a1bSJed Brown ym = info->ym; 1105c4762a1bSJed Brown zm = info->xm; 1106c4762a1bSJed Brown hx = thi->Lx / info->mz; 1107c4762a1bSJed Brown hy = thi->Ly / info->my; 1108c4762a1bSJed Brown 1109c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 1110c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 1111c4762a1bSJed Brown PrmNode pn[4],dpn[4][2]; 1112c4762a1bSJed Brown QuadExtract(prm,i,j,pn); 11139566063dSJacob Faibussowitsch PetscCall(QuadComputeGrad4(QuadQDeriv,hx,hy,pn,dpn)); 1114c4762a1bSJed Brown for (k=0; k<zm-1; k++) { 1115c4762a1bSJed Brown Node n[8]; 1116c4762a1bSJed Brown PetscReal zn[8],etabase = 0; 1117c4762a1bSJed Brown PetscScalar Ke[8*NODE_SIZE][8*NODE_SIZE],Kcpl[8*NODE_SIZE][4*PRMNODE_SIZE]; 1118c4762a1bSJed Brown PetscInt ls = 0; 1119c4762a1bSJed Brown 1120c4762a1bSJed Brown PrmHexGetZ(pn,k,zm,zn); 1121c4762a1bSJed Brown HexExtract(x,i,j,k,n); 11229566063dSJacob Faibussowitsch PetscCall(PetscMemzero(Ke,sizeof(Ke))); 11239566063dSJacob Faibussowitsch PetscCall(PetscMemzero(Kcpl,sizeof(Kcpl))); 1124c4762a1bSJed Brown if (thi->no_slip && k == 0) { 1125c4762a1bSJed Brown for (l=0; l<4; l++) n[l].u = n[l].v = 0; 1126c4762a1bSJed Brown ls = 4; 1127c4762a1bSJed Brown } 1128c4762a1bSJed Brown for (q=0; q<8; q++) { 1129c4762a1bSJed Brown PetscReal dz[3],phi[8],dphi[8][3],jw,eta,deta; 1130c4762a1bSJed Brown PetscScalar du[3],dv[3],u,v; 1131c4762a1bSJed Brown HexGrad(HexQDeriv[q],zn,dz); 1132c4762a1bSJed Brown HexComputeGeometry(q,hx,hy,dz,phi,dphi,&jw); 1133c4762a1bSJed Brown PointwiseNonlinearity(thi,n,phi,dphi,&u,&v,du,dv,&eta,&deta); 1134c4762a1bSJed Brown jw /= thi->rhog; /* residuals are scaled by this factor */ 1135c4762a1bSJed Brown if (q == 0) etabase = eta; 1136c4762a1bSJed Brown for (l=ls; l<8; l++) { /* test functions */ 1137c4762a1bSJed Brown const PetscReal pp=phi[l],*restrict dp = dphi[l]; 1138c4762a1bSJed Brown for (ll=ls; ll<8; ll++) { /* trial functions */ 1139c4762a1bSJed Brown const PetscReal *restrict dpl = dphi[ll]; 1140c4762a1bSJed Brown PetscScalar dgdu,dgdv; 1141c4762a1bSJed 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]; 1142c4762a1bSJed 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]; 1143c4762a1bSJed Brown /* Picard part */ 1144c4762a1bSJed 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]; 1145c4762a1bSJed Brown Ke[l*2+0][ll*2+1] += dp[0]*jw*eta*2.*dpl[1] + dp[1]*jw*eta*dpl[0]; 1146c4762a1bSJed Brown Ke[l*2+1][ll*2+0] += dp[1]*jw*eta*2.*dpl[0] + dp[0]*jw*eta*dpl[1]; 1147c4762a1bSJed 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]; 1148c4762a1bSJed Brown /* extra Newton terms */ 1149c4762a1bSJed 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]; 1150c4762a1bSJed 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]; 1151c4762a1bSJed 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]; 1152c4762a1bSJed 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]; 1153c4762a1bSJed Brown /* inertial part */ 1154c4762a1bSJed Brown Ke[l*2+0][ll*2+0] += pp*jw*thi->inertia*pp; 1155c4762a1bSJed Brown Ke[l*2+1][ll*2+1] += pp*jw*thi->inertia*pp; 1156c4762a1bSJed Brown } 1157c4762a1bSJed Brown for (ll=0; ll<4; ll++) { /* Trial functions for surface/bed */ 1158c4762a1bSJed Brown const PetscReal dpl[] = {QuadQDeriv[q%4][ll][0]/hx, QuadQDeriv[q%4][ll][1]/hy}; /* surface = h + b */ 1159c4762a1bSJed Brown Kcpl[FieldIndex(Node,l,u)][FieldIndex(PrmNode,ll,h)] += pp*jw*thi->rhog*dpl[0]; 1160c4762a1bSJed Brown Kcpl[FieldIndex(Node,l,u)][FieldIndex(PrmNode,ll,b)] += pp*jw*thi->rhog*dpl[0]; 1161c4762a1bSJed Brown Kcpl[FieldIndex(Node,l,v)][FieldIndex(PrmNode,ll,h)] += pp*jw*thi->rhog*dpl[1]; 1162c4762a1bSJed Brown Kcpl[FieldIndex(Node,l,v)][FieldIndex(PrmNode,ll,b)] += pp*jw*thi->rhog*dpl[1]; 1163c4762a1bSJed Brown } 1164c4762a1bSJed Brown } 1165c4762a1bSJed Brown } 1166c4762a1bSJed Brown if (k == 0) { /* on a bottom face */ 1167c4762a1bSJed Brown if (thi->no_slip) { 1168c4762a1bSJed Brown const PetscReal hz = PetscRealPart(pn[0].h)/(zm-1); 1169c4762a1bSJed 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); 1170c4762a1bSJed Brown Ke[0][0] = thi->dirichlet_scale*diagu; 1171c4762a1bSJed Brown Ke[0][1] = 0; 1172c4762a1bSJed Brown Ke[1][0] = 0; 1173c4762a1bSJed Brown Ke[1][1] = thi->dirichlet_scale*diagv; 1174c4762a1bSJed Brown } else { 1175c4762a1bSJed 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) */ 1176c4762a1bSJed Brown const PetscReal jw = 0.25*hx*hy,*phi = QuadQInterp[q]; 1177c4762a1bSJed Brown PetscScalar u =0,v=0,rbeta2=0; 1178c4762a1bSJed Brown PetscReal beta2,dbeta2; 1179c4762a1bSJed Brown for (l=0; l<4; l++) { 1180c4762a1bSJed Brown u += phi[l]*n[l].u; 1181c4762a1bSJed Brown v += phi[l]*n[l].v; 1182c4762a1bSJed Brown rbeta2 += phi[l]*pn[l].beta2; 1183c4762a1bSJed Brown } 1184c4762a1bSJed Brown THIFriction(thi,PetscRealPart(rbeta2),PetscRealPart(u*u+v*v)/2,&beta2,&dbeta2); 1185c4762a1bSJed Brown for (l=0; l<4; l++) { 1186c4762a1bSJed Brown const PetscReal pp = phi[l]; 1187c4762a1bSJed Brown for (ll=0; ll<4; ll++) { 1188c4762a1bSJed Brown const PetscReal ppl = phi[ll]; 1189c4762a1bSJed Brown Ke[l*2+0][ll*2+0] += pp*jw*beta2*ppl + pp*jw*dbeta2*u*u*ppl; 1190c4762a1bSJed Brown Ke[l*2+0][ll*2+1] += pp*jw*dbeta2*u*v*ppl; 1191c4762a1bSJed Brown Ke[l*2+1][ll*2+0] += pp*jw*dbeta2*v*u*ppl; 1192c4762a1bSJed Brown Ke[l*2+1][ll*2+1] += pp*jw*beta2*ppl + pp*jw*dbeta2*v*v*ppl; 1193c4762a1bSJed Brown } 1194c4762a1bSJed Brown } 1195c4762a1bSJed Brown } 1196c4762a1bSJed Brown } 1197c4762a1bSJed Brown } 1198c4762a1bSJed Brown { 1199c4762a1bSJed Brown const PetscInt rc3blocked[8] = { 1200c4762a1bSJed Brown DMDALocalIndex3D(info,i+0,j+0,k+0), 1201c4762a1bSJed Brown DMDALocalIndex3D(info,i+1,j+0,k+0), 1202c4762a1bSJed Brown DMDALocalIndex3D(info,i+1,j+1,k+0), 1203c4762a1bSJed Brown DMDALocalIndex3D(info,i+0,j+1,k+0), 1204c4762a1bSJed Brown DMDALocalIndex3D(info,i+0,j+0,k+1), 1205c4762a1bSJed Brown DMDALocalIndex3D(info,i+1,j+0,k+1), 1206c4762a1bSJed Brown DMDALocalIndex3D(info,i+1,j+1,k+1), 1207c4762a1bSJed Brown DMDALocalIndex3D(info,i+0,j+1,k+1) 1208c4762a1bSJed Brown },col2blocked[PRMNODE_SIZE*4] = { 1209c4762a1bSJed Brown DMDALocalIndex2D(info,i+0,j+0), 1210c4762a1bSJed Brown DMDALocalIndex2D(info,i+1,j+0), 1211c4762a1bSJed Brown DMDALocalIndex2D(info,i+1,j+1), 1212c4762a1bSJed Brown DMDALocalIndex2D(info,i+0,j+1) 1213c4762a1bSJed Brown }; 1214c4762a1bSJed Brown #if !defined COMPUTE_LOWER_TRIANGULAR /* fill in lower-triangular part, this is really cheap compared to computing the entries */ 1215c4762a1bSJed Brown for (l=0; l<8; l++) { 1216c4762a1bSJed Brown for (ll=l+1; ll<8; ll++) { 1217c4762a1bSJed Brown Ke[ll*2+0][l*2+0] = Ke[l*2+0][ll*2+0]; 1218c4762a1bSJed Brown Ke[ll*2+1][l*2+0] = Ke[l*2+0][ll*2+1]; 1219c4762a1bSJed Brown Ke[ll*2+0][l*2+1] = Ke[l*2+1][ll*2+0]; 1220c4762a1bSJed Brown Ke[ll*2+1][l*2+1] = Ke[l*2+1][ll*2+1]; 1221c4762a1bSJed Brown } 1222c4762a1bSJed Brown } 1223c4762a1bSJed Brown #endif 12249566063dSJacob Faibussowitsch PetscCall(MatSetValuesBlockedLocal(B,8,rc3blocked,8,rc3blocked,&Ke[0][0],ADD_VALUES)); /* velocity-velocity coupling can use blocked insertion */ 1225c4762a1bSJed Brown { /* The off-diagonal part cannot (yet) */ 1226c4762a1bSJed Brown PetscInt row3scalar[NODE_SIZE*8],col2scalar[PRMNODE_SIZE*4]; 1227c4762a1bSJed Brown for (l=0; l<8; l++) for (ll=0; ll<NODE_SIZE; ll++) row3scalar[l*NODE_SIZE+ll] = rc3blocked[l]*NODE_SIZE+ll; 1228c4762a1bSJed Brown for (l=0; l<4; l++) for (ll=0; ll<PRMNODE_SIZE; ll++) col2scalar[l*PRMNODE_SIZE+ll] = col2blocked[l]*PRMNODE_SIZE+ll; 12299566063dSJacob Faibussowitsch PetscCall(MatSetValuesLocal(Bcpl,8*NODE_SIZE,row3scalar,4*PRMNODE_SIZE,col2scalar,&Kcpl[0][0],ADD_VALUES)); 1230c4762a1bSJed Brown } 1231c4762a1bSJed Brown } 1232c4762a1bSJed Brown } 1233c4762a1bSJed Brown } 1234c4762a1bSJed Brown } 1235c4762a1bSJed Brown PetscFunctionReturn(0); 1236c4762a1bSJed Brown } 1237c4762a1bSJed Brown 1238c4762a1bSJed Brown static PetscErrorCode THIJacobianLocal_2D(DMDALocalInfo *info,const Node ***x3,const PrmNode **x2,const PrmNode **xdot2,PetscReal a,Mat B22,Mat B21,THI thi) 1239c4762a1bSJed Brown { 1240c4762a1bSJed Brown PetscInt xs,ys,xm,ym,zm,i,j,k; 1241c4762a1bSJed Brown 1242c4762a1bSJed Brown PetscFunctionBeginUser; 1243c4762a1bSJed Brown xs = info->zs; 1244c4762a1bSJed Brown ys = info->ys; 1245c4762a1bSJed Brown xm = info->zm; 1246c4762a1bSJed Brown ym = info->ym; 1247c4762a1bSJed Brown zm = info->xm; 1248c4762a1bSJed Brown 12493c633725SBarry Smith PetscCheck(zm <= 1024,((PetscObject)info->da)->comm,PETSC_ERR_SUP,"Need to allocate more space"); 1250c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 1251c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 1252c4762a1bSJed Brown { /* Self-coupling */ 1253c4762a1bSJed Brown const PetscInt row[] = {DMDALocalIndex2D(info,i,j)}; 1254c4762a1bSJed Brown const PetscInt col[] = {DMDALocalIndex2D(info,i,j)}; 1255c4762a1bSJed Brown const PetscScalar vals[] = { 1256c4762a1bSJed Brown a,0,0, 1257c4762a1bSJed Brown 0,a,0, 1258c4762a1bSJed Brown 0,0,a 1259c4762a1bSJed Brown }; 12609566063dSJacob Faibussowitsch PetscCall(MatSetValuesBlockedLocal(B22,1,row,1,col,vals,INSERT_VALUES)); 1261c4762a1bSJed Brown } 1262c4762a1bSJed Brown for (k=0; k<zm; k++) { /* Coupling to velocity problem */ 1263c4762a1bSJed Brown /* Use a cheaper quadrature than for residual evaluation, because it is much sparser */ 1264c4762a1bSJed Brown const PetscInt row[] = {FieldIndex(PrmNode,DMDALocalIndex2D(info,i,j),h)}; 1265c4762a1bSJed Brown const PetscInt cols[] = { 1266c4762a1bSJed Brown FieldIndex(Node,DMDALocalIndex3D(info,i-1,j,k),u), 1267c4762a1bSJed Brown FieldIndex(Node,DMDALocalIndex3D(info,i ,j,k),u), 1268c4762a1bSJed Brown FieldIndex(Node,DMDALocalIndex3D(info,i+1,j,k),u), 1269c4762a1bSJed Brown FieldIndex(Node,DMDALocalIndex3D(info,i,j-1,k),v), 1270c4762a1bSJed Brown FieldIndex(Node,DMDALocalIndex3D(info,i,j ,k),v), 1271c4762a1bSJed Brown FieldIndex(Node,DMDALocalIndex3D(info,i,j+1,k),v) 1272c4762a1bSJed Brown }; 1273c4762a1bSJed Brown const PetscScalar 1274c4762a1bSJed Brown w = (k && k<zm-1) ? 0.5 : 0.25, 1275c4762a1bSJed Brown hW = w*(x2[i-1][j ].h+x2[i ][j ].h)/(zm-1.), 1276c4762a1bSJed Brown hE = w*(x2[i ][j ].h+x2[i+1][j ].h)/(zm-1.), 1277c4762a1bSJed Brown hS = w*(x2[i ][j-1].h+x2[i ][j ].h)/(zm-1.), 1278c4762a1bSJed Brown hN = w*(x2[i ][j ].h+x2[i ][j+1].h)/(zm-1.); 1279c4762a1bSJed Brown PetscScalar *vals, 1280c4762a1bSJed Brown vals_upwind[] = {((PetscRealPart(x3[i][j][k].u) > 0) ? -hW : 0), 1281c4762a1bSJed Brown ((PetscRealPart(x3[i][j][k].u) > 0) ? +hE : -hW), 1282c4762a1bSJed Brown ((PetscRealPart(x3[i][j][k].u) > 0) ? 0 : +hE), 1283c4762a1bSJed Brown ((PetscRealPart(x3[i][j][k].v) > 0) ? -hS : 0), 1284c4762a1bSJed Brown ((PetscRealPart(x3[i][j][k].v) > 0) ? +hN : -hS), 1285c4762a1bSJed Brown ((PetscRealPart(x3[i][j][k].v) > 0) ? 0 : +hN)}, 1286c4762a1bSJed Brown vals_centered[] = {-0.5*hW, 0.5*(-hW+hE), 0.5*hE, 1287c4762a1bSJed Brown -0.5*hS, 0.5*(-hS+hN), 0.5*hN}; 1288c4762a1bSJed Brown vals = 1 ? vals_upwind : vals_centered; 1289c4762a1bSJed Brown if (k == 0) { 1290c4762a1bSJed Brown Node derate; 1291c4762a1bSJed Brown THIErosion(thi,&x3[i][j][0],NULL,&derate); 1292c4762a1bSJed Brown vals[1] -= derate.u; 1293c4762a1bSJed Brown vals[4] -= derate.v; 1294c4762a1bSJed Brown } 12959566063dSJacob Faibussowitsch PetscCall(MatSetValuesLocal(B21,1,row,6,cols,vals,INSERT_VALUES)); 1296c4762a1bSJed Brown } 1297c4762a1bSJed Brown } 1298c4762a1bSJed Brown } 1299c4762a1bSJed Brown PetscFunctionReturn(0); 1300c4762a1bSJed Brown } 1301c4762a1bSJed Brown 1302c4762a1bSJed Brown static PetscErrorCode THIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,void *ctx) 1303c4762a1bSJed Brown { 1304c4762a1bSJed Brown THI thi = (THI)ctx; 1305c4762a1bSJed Brown DM pack,da3,da2; 1306c4762a1bSJed Brown Vec X3,X2,Xdot2; 1307c4762a1bSJed Brown Mat B11,B12,B21,B22; 1308c4762a1bSJed Brown DMDALocalInfo info3; 1309c4762a1bSJed Brown IS *isloc; 1310c4762a1bSJed Brown const Node ***x3; 1311c4762a1bSJed Brown const PrmNode **x2,**xdot2; 1312c4762a1bSJed Brown 1313c4762a1bSJed Brown PetscFunctionBeginUser; 13149566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts,&pack)); 13159566063dSJacob Faibussowitsch PetscCall(DMCompositeGetEntries(pack,&da3,&da2)); 13169566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da3,&info3)); 13179566063dSJacob Faibussowitsch PetscCall(DMCompositeGetLocalVectors(pack,&X3,&X2)); 13189566063dSJacob Faibussowitsch PetscCall(DMCompositeGetLocalVectors(pack,NULL,&Xdot2)); 13199566063dSJacob Faibussowitsch PetscCall(DMCompositeScatter(pack,X,X3,X2)); 13209566063dSJacob Faibussowitsch PetscCall(THIFixGhosts(thi,da3,da2,X3,X2)); 13219566063dSJacob Faibussowitsch PetscCall(DMCompositeScatter(pack,Xdot,NULL,Xdot2)); 1322c4762a1bSJed Brown 13239566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(B)); 1324c4762a1bSJed Brown 13259566063dSJacob Faibussowitsch PetscCall(DMCompositeGetLocalISs(pack,&isloc)); 13269566063dSJacob Faibussowitsch PetscCall(MatGetLocalSubMatrix(B,isloc[0],isloc[0],&B11)); 13279566063dSJacob Faibussowitsch PetscCall(MatGetLocalSubMatrix(B,isloc[0],isloc[1],&B12)); 13289566063dSJacob Faibussowitsch PetscCall(MatGetLocalSubMatrix(B,isloc[1],isloc[0],&B21)); 13299566063dSJacob Faibussowitsch PetscCall(MatGetLocalSubMatrix(B,isloc[1],isloc[1],&B22)); 1330c4762a1bSJed Brown 13319566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da3,X3,&x3)); 13329566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da2,X2,&x2)); 13339566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da2,Xdot2,&xdot2)); 1334c4762a1bSJed Brown 13359566063dSJacob Faibussowitsch PetscCall(THIJacobianLocal_Momentum(&info3,x3,x2,B11,B12,thi)); 1336c4762a1bSJed Brown 1337c4762a1bSJed Brown /* Need to switch from ADD_VALUES to INSERT_VALUES */ 13389566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B,MAT_FLUSH_ASSEMBLY)); 13399566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B,MAT_FLUSH_ASSEMBLY)); 1340c4762a1bSJed Brown 13419566063dSJacob Faibussowitsch PetscCall(THIJacobianLocal_2D(&info3,x3,x2,xdot2,a,B22,B21,thi)); 1342c4762a1bSJed Brown 13439566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da3,X3,&x3)); 13449566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da2,X2,&x2)); 13459566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da2,Xdot2,&xdot2)); 1346c4762a1bSJed Brown 13479566063dSJacob Faibussowitsch PetscCall(MatRestoreLocalSubMatrix(B,isloc[0],isloc[0],&B11)); 13489566063dSJacob Faibussowitsch PetscCall(MatRestoreLocalSubMatrix(B,isloc[0],isloc[1],&B12)); 13499566063dSJacob Faibussowitsch PetscCall(MatRestoreLocalSubMatrix(B,isloc[1],isloc[0],&B21)); 13509566063dSJacob Faibussowitsch PetscCall(MatRestoreLocalSubMatrix(B,isloc[1],isloc[1],&B22)); 13519566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isloc[0])); 13529566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isloc[1])); 13539566063dSJacob Faibussowitsch PetscCall(PetscFree(isloc)); 1354c4762a1bSJed Brown 13559566063dSJacob Faibussowitsch PetscCall(DMCompositeRestoreLocalVectors(pack,&X3,&X2)); 13569566063dSJacob Faibussowitsch PetscCall(DMCompositeRestoreLocalVectors(pack,0,&Xdot2)); 1357c4762a1bSJed Brown 13589566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 13599566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 1360c4762a1bSJed Brown if (A != B) { 13619566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 13629566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 1363c4762a1bSJed Brown } 13649566063dSJacob Faibussowitsch if (thi->verbose) PetscCall(THIMatrixStatistics(thi,B,PETSC_VIEWER_STDOUT_WORLD)); 1365c4762a1bSJed Brown PetscFunctionReturn(0); 1366c4762a1bSJed Brown } 1367c4762a1bSJed Brown 1368c4762a1bSJed Brown /* VTK's XML formats are so brain-dead that they can't handle multiple grids in the same file. Since the communication 1369c4762a1bSJed Brown * can be shared between the two grids, we write two files at once, one for velocity (living on a 3D grid defined by 1370c4762a1bSJed Brown * h=thickness and b=bed) and another for all properties living on the 2D grid. 1371c4762a1bSJed Brown */ 1372c4762a1bSJed Brown static PetscErrorCode THIDAVecView_VTK_XML(THI thi,DM pack,Vec X,const char filename[],const char filename2[]) 1373c4762a1bSJed Brown { 1374c4762a1bSJed Brown const PetscInt dof = NODE_SIZE,dof2 = PRMNODE_SIZE; 1375c4762a1bSJed Brown Units units = thi->units; 1376c4762a1bSJed Brown MPI_Comm comm; 1377c4762a1bSJed Brown PetscViewer viewer3,viewer2; 1378c4762a1bSJed Brown PetscMPIInt rank,size,tag,nn,nmax,nn2,nmax2; 1379c4762a1bSJed Brown PetscInt mx,my,mz,r,range[6]; 1380c4762a1bSJed Brown PetscScalar *x,*x2; 1381c4762a1bSJed Brown DM da3,da2; 1382c4762a1bSJed Brown Vec X3,X2; 1383c4762a1bSJed Brown 1384c4762a1bSJed Brown PetscFunctionBeginUser; 13859566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)thi,&comm)); 13869566063dSJacob Faibussowitsch PetscCall(DMCompositeGetEntries(pack,&da3,&da2)); 13879566063dSJacob Faibussowitsch PetscCall(DMCompositeGetAccess(pack,X,&X3,&X2)); 13889566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da3,0, &mz,&my,&mx, 0,0,0, 0,0,0,0,0,0)); 13899566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm,&size)); 13909566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm,&rank)); 13919566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIOpen(comm,filename,&viewer3)); 13929566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIOpen(comm,filename2,&viewer2)); 13939566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer3,"<VTKFile type=\"StructuredGrid\" version=\"0.1\" byte_order=\"LittleEndian\">\n")); 13949566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer2,"<VTKFile type=\"StructuredGrid\" version=\"0.1\" byte_order=\"LittleEndian\">\n")); 139563a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer3," <StructuredGrid WholeExtent=\"%d %" PetscInt_FMT " %d %" PetscInt_FMT " %d %" PetscInt_FMT "\">\n",0,mz-1,0,my-1,0,mx-1)); 139663a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer2," <StructuredGrid WholeExtent=\"%d %d %d %" PetscInt_FMT " %d %" PetscInt_FMT "\">\n",0,0,0,my-1,0,mx-1)); 1397c4762a1bSJed Brown 13989566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da3,range,range+1,range+2,range+3,range+4,range+5)); 13999566063dSJacob Faibussowitsch PetscCall(PetscMPIIntCast(range[3]*range[4]*range[5]*dof,&nn)); 14009566063dSJacob Faibussowitsch PetscCallMPI(MPI_Reduce(&nn,&nmax,1,MPI_INT,MPI_MAX,0,comm)); 14019566063dSJacob Faibussowitsch PetscCall(PetscMPIIntCast(range[4]*range[5]*dof2,&nn2)); 14029566063dSJacob Faibussowitsch PetscCallMPI(MPI_Reduce(&nn2,&nmax2,1,MPI_INT,MPI_MAX,0,comm)); 1403c4762a1bSJed Brown tag = ((PetscObject)viewer3)->tag; 14049566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X3,(const PetscScalar**)&x)); 14059566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X2,(const PetscScalar**)&x2)); 1406dd400576SPatrick Sanan if (rank == 0) { 1407c4762a1bSJed Brown PetscScalar *array,*array2; 14089566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(nmax,&array,nmax2,&array2)); 1409c4762a1bSJed Brown for (r=0; r<size; r++) { 1410c4762a1bSJed Brown PetscInt i,j,k,f,xs,xm,ys,ym,zs,zm; 1411c4762a1bSJed Brown Node *y3; 1412c4762a1bSJed Brown PetscScalar (*y2)[PRMNODE_SIZE]; 1413c4762a1bSJed Brown MPI_Status status; 1414c4762a1bSJed Brown 1415c4762a1bSJed Brown if (r) { 14169566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(range,6,MPIU_INT,r,tag,comm,MPI_STATUS_IGNORE)); 1417c4762a1bSJed Brown } 1418c4762a1bSJed Brown zs = range[0];ys = range[1];xs = range[2];zm = range[3];ym = range[4];xm = range[5]; 14193c633725SBarry Smith PetscCheck(xm*ym*zm*dof <= nmax,PETSC_COMM_SELF,PETSC_ERR_PLIB,"should not happen"); 1420c4762a1bSJed Brown if (r) { 14219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(array,nmax,MPIU_SCALAR,r,tag,comm,&status)); 14229566063dSJacob Faibussowitsch PetscCallMPI(MPI_Get_count(&status,MPIU_SCALAR,&nn)); 14233c633725SBarry Smith PetscCheck(nn == xm*ym*zm*dof,PETSC_COMM_SELF,PETSC_ERR_PLIB,"corrupt da3 send"); 1424c4762a1bSJed Brown y3 = (Node*)array; 14259566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(array2,nmax2,MPIU_SCALAR,r,tag,comm,&status)); 14269566063dSJacob Faibussowitsch PetscCallMPI(MPI_Get_count(&status,MPIU_SCALAR,&nn2)); 14273c633725SBarry Smith PetscCheck(nn2 == xm*ym*dof2,PETSC_COMM_SELF,PETSC_ERR_PLIB,"corrupt da2 send"); 1428c4762a1bSJed Brown y2 = (PetscScalar(*)[PRMNODE_SIZE])array2; 1429c4762a1bSJed Brown } else { 1430c4762a1bSJed Brown y3 = (Node*)x; 1431c4762a1bSJed Brown y2 = (PetscScalar(*)[PRMNODE_SIZE])x2; 1432c4762a1bSJed Brown } 143363a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer3," <Piece Extent=\"%" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\">\n",zs,zs+zm-1,ys,ys+ym-1,xs,xs+xm-1)); 143463a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer2," <Piece Extent=\"%d %d %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\">\n",0,0,ys,ys+ym-1,xs,xs+xm-1)); 1435c4762a1bSJed Brown 14369566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer3," <Points>\n")); 14379566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer2," <Points>\n")); 14389566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer3," <DataArray type=\"Float32\" NumberOfComponents=\"3\" format=\"ascii\">\n")); 14399566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer2," <DataArray type=\"Float32\" NumberOfComponents=\"3\" format=\"ascii\">\n")); 1440c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 1441c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 1442c4762a1bSJed Brown PetscReal 1443c4762a1bSJed Brown xx = thi->Lx*i/mx, 1444c4762a1bSJed Brown yy = thi->Ly*j/my, 1445c4762a1bSJed Brown b = PetscRealPart(y2[i*ym+j][FieldOffset(PrmNode,b)]), 1446c4762a1bSJed Brown h = PetscRealPart(y2[i*ym+j][FieldOffset(PrmNode,h)]); 1447c4762a1bSJed Brown for (k=zs; k<zs+zm; k++) { 1448c4762a1bSJed Brown PetscReal zz = b + h*k/(mz-1.); 144963a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer3,"%f %f %f\n",(double)xx,(double)yy,(double)zz)); 1450c4762a1bSJed Brown } 145163a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer2,"%f %f %f\n",(double)xx,(double)yy,(double)0.0)); 1452c4762a1bSJed Brown } 1453c4762a1bSJed Brown } 14549566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer3," </DataArray>\n")); 14559566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer2," </DataArray>\n")); 14569566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer3," </Points>\n")); 14579566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer2," </Points>\n")); 1458c4762a1bSJed Brown 1459c4762a1bSJed Brown { /* Velocity and rank (3D) */ 14609566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer3," <PointData>\n")); 14619566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer3," <DataArray type=\"Float32\" Name=\"velocity\" NumberOfComponents=\"3\" format=\"ascii\">\n")); 1462c4762a1bSJed Brown for (i=0; i<nn/dof; i++) { 146363a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer3,"%f %f %f\n",(double)(PetscRealPart(y3[i].u)*units->year/units->meter),(double)(PetscRealPart(y3[i].v)*units->year/units->meter),0.0)); 1464c4762a1bSJed Brown } 14659566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer3," </DataArray>\n")); 1466c4762a1bSJed Brown 14679566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer3," <DataArray type=\"Int32\" Name=\"rank\" NumberOfComponents=\"1\" format=\"ascii\">\n")); 1468c4762a1bSJed Brown for (i=0; i<nn; i+=dof) { 146963a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer3,"%" PetscInt_FMT "\n",r)); 1470c4762a1bSJed Brown } 14719566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer3," </DataArray>\n")); 14729566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer3," </PointData>\n")); 1473c4762a1bSJed Brown } 1474c4762a1bSJed Brown 1475c4762a1bSJed Brown { /* 2D */ 14769566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer2," <PointData>\n")); 1477c4762a1bSJed Brown for (f=0; f<PRMNODE_SIZE; f++) { 1478c4762a1bSJed Brown const char *fieldname; 14799566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldName(da2,f,&fieldname)); 14809566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer2," <DataArray type=\"Float32\" Name=\"%s\" format=\"ascii\">\n",fieldname)); 1481c4762a1bSJed Brown for (i=0; i<nn2/PRMNODE_SIZE; i++) { 148263a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer2,"%g\n",(double)y2[i][f])); 1483c4762a1bSJed Brown } 14849566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer2," </DataArray>\n")); 1485c4762a1bSJed Brown } 14869566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer2," </PointData>\n")); 1487c4762a1bSJed Brown } 1488c4762a1bSJed Brown 14899566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer3," </Piece>\n")); 14909566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer2," </Piece>\n")); 1491c4762a1bSJed Brown } 14929566063dSJacob Faibussowitsch PetscCall(PetscFree2(array,array2)); 1493c4762a1bSJed Brown } else { 14949566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send(range,6,MPIU_INT,0,tag,comm)); 14959566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send(x,nn,MPIU_SCALAR,0,tag,comm)); 14969566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send(x2,nn2,MPIU_SCALAR,0,tag,comm)); 1497c4762a1bSJed Brown } 14989566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X3,(const PetscScalar**)&x)); 14999566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X2,(const PetscScalar**)&x2)); 15009566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer3," </StructuredGrid>\n")); 15019566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer2," </StructuredGrid>\n")); 1502c4762a1bSJed Brown 15039566063dSJacob Faibussowitsch PetscCall(DMCompositeRestoreAccess(pack,X,&X3,&X2)); 15049566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer3,"</VTKFile>\n")); 15059566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer2,"</VTKFile>\n")); 15069566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer3)); 15079566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer2)); 1508c4762a1bSJed Brown PetscFunctionReturn(0); 1509c4762a1bSJed Brown } 1510c4762a1bSJed Brown 1511c4762a1bSJed Brown static PetscErrorCode THITSMonitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx) 1512c4762a1bSJed Brown { 1513c4762a1bSJed Brown THI thi = (THI)ctx; 1514c4762a1bSJed Brown DM pack; 1515c4762a1bSJed Brown char filename3[PETSC_MAX_PATH_LEN],filename2[PETSC_MAX_PATH_LEN]; 1516c4762a1bSJed Brown 1517c4762a1bSJed Brown PetscFunctionBeginUser; 1518c4762a1bSJed Brown if (step < 0) PetscFunctionReturn(0); /* negative one is used to indicate an interpolated solution */ 151963a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts),"%3" PetscInt_FMT ": t=%g\n",step,(double)t)); 1520c4762a1bSJed Brown if (thi->monitor_interval && step % thi->monitor_interval) PetscFunctionReturn(0); 15219566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts,&pack)); 152263a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(filename3,sizeof(filename3),"%s-3d-%03" PetscInt_FMT ".vts",thi->monitor_basename,step)); 152363a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(filename2,sizeof(filename2),"%s-2d-%03" PetscInt_FMT ".vts",thi->monitor_basename,step)); 15249566063dSJacob Faibussowitsch PetscCall(THIDAVecView_VTK_XML(thi,pack,X,filename3,filename2)); 1525c4762a1bSJed Brown PetscFunctionReturn(0); 1526c4762a1bSJed Brown } 1527c4762a1bSJed Brown 1528c4762a1bSJed Brown static PetscErrorCode THICreateDM3d(THI thi,DM *dm3d) 1529c4762a1bSJed Brown { 1530c4762a1bSJed Brown MPI_Comm comm; 1531c4762a1bSJed Brown PetscInt M = 3,N = 3,P = 2; 1532c4762a1bSJed Brown DM da; 1533c4762a1bSJed Brown 1534c4762a1bSJed Brown PetscFunctionBeginUser; 15359566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)thi,&comm)); 1536d0609cedSBarry Smith PetscOptionsBegin(comm,NULL,"Grid resolution options",""); 1537c4762a1bSJed Brown { 15389566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-M","Number of elements in x-direction on coarse level","",M,&M,NULL)); 1539c4762a1bSJed Brown N = M; 15409566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-N","Number of elements in y-direction on coarse level (if different from M)","",N,&N,NULL)); 15419566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-P","Number of elements in z-direction on coarse level","",P,&P,NULL)); 1542c4762a1bSJed Brown } 1543d0609cedSBarry Smith PetscOptionsEnd(); 15449566063dSJacob Faibussowitsch PetscCall(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)); 15459566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 15469566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 15479566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da,0,"x-velocity")); 15489566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da,1,"y-velocity")); 1549c4762a1bSJed Brown *dm3d = da; 1550c4762a1bSJed Brown PetscFunctionReturn(0); 1551c4762a1bSJed Brown } 1552c4762a1bSJed Brown 1553c4762a1bSJed Brown int main(int argc,char *argv[]) 1554c4762a1bSJed Brown { 1555c4762a1bSJed Brown MPI_Comm comm; 1556c4762a1bSJed Brown DM pack,da3,da2; 1557c4762a1bSJed Brown TS ts; 1558c4762a1bSJed Brown THI thi; 1559c4762a1bSJed Brown Vec X; 1560c4762a1bSJed Brown Mat B; 1561c4762a1bSJed Brown PetscInt i,steps; 1562c4762a1bSJed Brown PetscReal ftime; 1563c4762a1bSJed Brown 1564*327415f7SBarry Smith PetscFunctionBeginUser; 15659566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,0,help)); 1566c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 1567c4762a1bSJed Brown 15689566063dSJacob Faibussowitsch PetscCall(THICreate(comm,&thi)); 15699566063dSJacob Faibussowitsch PetscCall(THICreateDM3d(thi,&da3)); 1570c4762a1bSJed Brown { 1571c4762a1bSJed Brown PetscInt Mx,My,mx,my,s; 1572c4762a1bSJed Brown DMDAStencilType st; 15739566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da3,0, 0,&My,&Mx, 0,&my,&mx, 0,&s,0,0,0,&st)); 15749566063dSJacob Faibussowitsch PetscCall(DMDACreate2d(PetscObjectComm((PetscObject)thi),DM_BOUNDARY_PERIODIC,DM_BOUNDARY_PERIODIC,st,My,Mx,my,mx,sizeof(PrmNode)/sizeof(PetscScalar),s,0,0,&da2)); 15759566063dSJacob Faibussowitsch PetscCall(DMSetUp(da2)); 1576c4762a1bSJed Brown } 1577c4762a1bSJed Brown 15789566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)da3,"3D_Velocity")); 15799566063dSJacob Faibussowitsch PetscCall(DMSetOptionsPrefix(da3,"f3d_")); 15809566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da3,0,"u")); 15819566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da3,1,"v")); 15829566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)da2,"2D_Fields")); 15839566063dSJacob Faibussowitsch PetscCall(DMSetOptionsPrefix(da2,"f2d_")); 15849566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da2,0,"b")); 15859566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da2,1,"h")); 15869566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da2,2,"beta2")); 15879566063dSJacob Faibussowitsch PetscCall(DMCompositeCreate(comm,&pack)); 15889566063dSJacob Faibussowitsch PetscCall(DMCompositeAddDM(pack,da3)); 15899566063dSJacob Faibussowitsch PetscCall(DMCompositeAddDM(pack,da2)); 15909566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da3)); 15919566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da2)); 15929566063dSJacob Faibussowitsch PetscCall(DMSetUp(pack)); 15939566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(pack,&B)); 15949566063dSJacob Faibussowitsch PetscCall(MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE)); 15959566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(B,"thi_")); 1596c4762a1bSJed Brown 1597c4762a1bSJed Brown for (i=0; i<thi->nlevels; i++) { 1598c4762a1bSJed Brown PetscReal Lx = thi->Lx / thi->units->meter,Ly = thi->Ly / thi->units->meter,Lz = thi->Lz / thi->units->meter; 1599c4762a1bSJed Brown PetscInt Mx,My,Mz; 16009566063dSJacob Faibussowitsch PetscCall(DMCompositeGetEntries(pack,&da3,&da2)); 16019566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da3,0, &Mz,&My,&Mx, 0,0,0, 0,0,0,0,0,0)); 160263a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)thi),"Level %" PetscInt_FMT " 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))); 1603c4762a1bSJed Brown } 1604c4762a1bSJed Brown 16059566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(pack,&X)); 16069566063dSJacob Faibussowitsch PetscCall(THIInitial(thi,pack,X)); 1607c4762a1bSJed Brown 16089566063dSJacob Faibussowitsch PetscCall(TSCreate(comm,&ts)); 16099566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts,pack)); 16109566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts,TS_NONLINEAR)); 16119566063dSJacob Faibussowitsch PetscCall(TSMonitorSet(ts,THITSMonitor,thi,NULL)); 16129566063dSJacob Faibussowitsch PetscCall(TSSetType(ts,TSTHETA)); 16139566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts,NULL,THIFunction,thi)); 16149566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(ts,B,B,THIJacobian,thi)); 16159566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts,10.0)); 16169566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER)); 16179566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts,X)); 16189566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts,1e-3)); 16199566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 1620c4762a1bSJed Brown 16219566063dSJacob Faibussowitsch PetscCall(TSSolve(ts,X)); 16229566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts,&ftime)); 16239566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts,&steps)); 162463a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Steps %" PetscInt_FMT " final time %g\n",steps,(double)ftime)); 1625c4762a1bSJed Brown 16269566063dSJacob Faibussowitsch if (0) PetscCall(THISolveStatistics(thi,ts,0,"Full")); 1627c4762a1bSJed Brown 1628c4762a1bSJed Brown { 1629c4762a1bSJed Brown PetscBool flg; 1630c4762a1bSJed Brown char filename[PETSC_MAX_PATH_LEN] = ""; 16319566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL,NULL,"-o",filename,sizeof(filename),&flg)); 16321baa6e33SBarry Smith if (flg) PetscCall(THIDAVecView_VTK_XML(thi,pack,X,filename,NULL)); 1633c4762a1bSJed Brown } 1634c4762a1bSJed Brown 16359566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X)); 16369566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 16379566063dSJacob Faibussowitsch PetscCall(DMDestroy(&pack)); 16389566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 16399566063dSJacob Faibussowitsch PetscCall(THIDestroy(&thi)); 16409566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 1641b122ec5aSJacob Faibussowitsch return 0; 1642c4762a1bSJed Brown } 1643