1af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 337c2070cSMatthew G. Knepley #include <petsc/private/snesimpl.h> 4924a1b8fSMatthew G. Knepley #include <petscds.h> 56dbbd306SMatthew G. Knepley #include <petscfv.h> 66dbbd306SMatthew G. Knepley 76da023fcSToby Isaac static PetscErrorCode DMTSConvertPlex(DM dm, DM *plex, PetscBool copy) 86da023fcSToby Isaac { 96da023fcSToby Isaac PetscBool isPlex; 106da023fcSToby Isaac PetscErrorCode ierr; 116da023fcSToby Isaac 126da023fcSToby Isaac PetscFunctionBegin; 136da023fcSToby Isaac ierr = PetscObjectTypeCompare((PetscObject) dm, DMPLEX, &isPlex);CHKERRQ(ierr); 146da023fcSToby Isaac if (isPlex) { 156da023fcSToby Isaac *plex = dm; 166da023fcSToby Isaac ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr); 17f7148743SMatthew G. Knepley } else { 18f7148743SMatthew G. Knepley ierr = PetscObjectQuery((PetscObject) dm, "dm_plex", (PetscObject *) plex);CHKERRQ(ierr); 19f7148743SMatthew G. Knepley if (!*plex) { 206da023fcSToby Isaac ierr = DMConvert(dm,DMPLEX,plex);CHKERRQ(ierr); 21f7148743SMatthew G. Knepley ierr = PetscObjectCompose((PetscObject) dm, "dm_plex", (PetscObject) *plex);CHKERRQ(ierr); 226da023fcSToby Isaac if (copy) { 236da023fcSToby Isaac PetscInt i; 246da023fcSToby Isaac PetscObject obj; 256da023fcSToby Isaac const char *comps[3] = {"A","dmAux","dmCh"}; 266da023fcSToby Isaac 276da023fcSToby Isaac ierr = DMCopyDMTS(dm, *plex);CHKERRQ(ierr); 286da023fcSToby Isaac ierr = DMCopyDMSNES(dm, *plex);CHKERRQ(ierr); 296da023fcSToby Isaac for (i = 0; i < 3; i++) { 306da023fcSToby Isaac ierr = PetscObjectQuery((PetscObject) dm, comps[i], &obj);CHKERRQ(ierr); 316da023fcSToby Isaac ierr = PetscObjectCompose((PetscObject) *plex, comps[i], obj);CHKERRQ(ierr); 326da023fcSToby Isaac } 336da023fcSToby Isaac } 34f7148743SMatthew G. Knepley } else { 35f7148743SMatthew G. Knepley ierr = PetscObjectReference((PetscObject) *plex);CHKERRQ(ierr); 36f7148743SMatthew G. Knepley } 376da023fcSToby Isaac } 386da023fcSToby Isaac PetscFunctionReturn(0); 396da023fcSToby Isaac } 406da023fcSToby Isaac 4108449791SMatthew G. Knepley /*@ 4208449791SMatthew G. Knepley DMPlexTSComputeRHSFunctionFVM - Form the local forcing F from the local input X using pointwise functions specified by the user 4308449791SMatthew G. Knepley 4408449791SMatthew G. Knepley Input Parameters: 4508449791SMatthew G. Knepley + dm - The mesh 4608449791SMatthew G. Knepley . t - The time 4708449791SMatthew G. Knepley . locX - Local solution 4808449791SMatthew G. Knepley - user - The user context 4908449791SMatthew G. Knepley 5008449791SMatthew G. Knepley Output Parameter: 513b16df42SMatthew G. Knepley . F - Global output vector 5208449791SMatthew G. Knepley 5308449791SMatthew G. Knepley Level: developer 5408449791SMatthew G. Knepley 5508449791SMatthew G. Knepley .seealso: DMPlexComputeJacobianActionFEM() 5608449791SMatthew G. Knepley @*/ 573b16df42SMatthew G. Knepley PetscErrorCode DMPlexTSComputeRHSFunctionFVM(DM dm, PetscReal time, Vec locX, Vec F, void *user) 58254c1ad2SMatthew G. Knepley { 593b16df42SMatthew G. Knepley Vec locF; 604a3e9fdbSToby Isaac IS cellIS; 616da023fcSToby Isaac DM plex; 624a3e9fdbSToby Isaac PetscInt depth; 63254c1ad2SMatthew G. Knepley PetscErrorCode ierr; 64254c1ad2SMatthew G. Knepley 65254c1ad2SMatthew G. Knepley PetscFunctionBegin; 666da023fcSToby Isaac ierr = DMTSConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr); 674a3e9fdbSToby Isaac ierr = DMPlexGetDepth(plex, &depth);CHKERRQ(ierr); 68aeadca18SToby Isaac ierr = DMGetStratumIS(plex, "dim", depth, &cellIS);CHKERRQ(ierr); 694a3e9fdbSToby Isaac if (!cellIS) { 704a3e9fdbSToby Isaac ierr = DMGetStratumIS(plex, "depth", depth, &cellIS);CHKERRQ(ierr); 714a3e9fdbSToby Isaac } 726da023fcSToby Isaac ierr = DMGetLocalVector(plex, &locF);CHKERRQ(ierr); 733b16df42SMatthew G. Knepley ierr = VecZeroEntries(locF);CHKERRQ(ierr); 744a3e9fdbSToby Isaac ierr = DMPlexComputeResidual_Internal(plex, cellIS, time, locX, NULL, time, locF, user);CHKERRQ(ierr); 759bda831aSToby Isaac ierr = DMLocalToGlobalBegin(plex, locF, ADD_VALUES, F);CHKERRQ(ierr); 769bda831aSToby Isaac ierr = DMLocalToGlobalEnd(plex, locF, ADD_VALUES, F);CHKERRQ(ierr); 776da023fcSToby Isaac ierr = DMRestoreLocalVector(plex, &locF);CHKERRQ(ierr); 784a3e9fdbSToby Isaac ierr = ISDestroy(&cellIS);CHKERRQ(ierr); 79b2338dbcSToby Isaac ierr = DMDestroy(&plex);CHKERRQ(ierr); 80254c1ad2SMatthew G. Knepley PetscFunctionReturn(0); 81254c1ad2SMatthew G. Knepley } 82254c1ad2SMatthew G. Knepley 83c5d70e09SMatthew G. Knepley /*@ 84ef68eab9SMatthew G. Knepley DMPlexTSComputeBoundary - Insert the essential boundary values for the local input X and/or its time derivative X_t using pointwise functions specified by the user 85c5d70e09SMatthew G. Knepley 86c5d70e09SMatthew G. Knepley Input Parameters: 87c5d70e09SMatthew G. Knepley + dm - The mesh 88a40652d4SToby Isaac . t - The time 89a40652d4SToby Isaac . locX - Local solution 90a40652d4SToby Isaac . locX_t - Local solution time derivative, or NULL 91c5d70e09SMatthew G. Knepley - user - The user context 92c5d70e09SMatthew G. Knepley 93c5d70e09SMatthew G. Knepley Level: developer 94c5d70e09SMatthew G. Knepley 95a40652d4SToby Isaac .seealso: DMPlexComputeJacobianActionFEM() 96c5d70e09SMatthew G. Knepley @*/ 97ef68eab9SMatthew G. Knepley PetscErrorCode DMPlexTSComputeBoundary(DM dm, PetscReal time, Vec locX, Vec locX_t, void *user) 98c5d70e09SMatthew G. Knepley { 99c5d70e09SMatthew G. Knepley DM plex; 100ef68eab9SMatthew G. Knepley Vec faceGeometryFVM = NULL; 101ef68eab9SMatthew G. Knepley PetscInt Nf, f; 102c5d70e09SMatthew G. Knepley PetscErrorCode ierr; 103c5d70e09SMatthew G. Knepley 104c5d70e09SMatthew G. Knepley PetscFunctionBegin; 105c5d70e09SMatthew G. Knepley ierr = DMTSConvertPlex(dm, &plex, PETSC_TRUE);CHKERRQ(ierr); 106ef68eab9SMatthew G. Knepley ierr = DMGetNumFields(plex, &Nf);CHKERRQ(ierr); 1079586001cSMatthew G. Knepley if (!locX_t) { 1089586001cSMatthew G. Knepley /* This is the RHS part */ 109ef68eab9SMatthew G. Knepley for (f = 0; f < Nf; f++) { 110ef68eab9SMatthew G. Knepley PetscObject obj; 111ef68eab9SMatthew G. Knepley PetscClassId id; 112ef68eab9SMatthew G. Knepley 11344a7f3ddSMatthew G. Knepley ierr = DMGetField(plex, f, NULL, &obj);CHKERRQ(ierr); 114ef68eab9SMatthew G. Knepley ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 115ef68eab9SMatthew G. Knepley if (id == PETSCFV_CLASSID) { 1163e9753d6SMatthew G. Knepley ierr = DMPlexGetGeometryFVM(plex, &faceGeometryFVM, NULL, NULL);CHKERRQ(ierr); 117ef68eab9SMatthew G. Knepley break; 118ef68eab9SMatthew G. Knepley } 119ef68eab9SMatthew G. Knepley } 1209586001cSMatthew G. Knepley } 1219586001cSMatthew G. Knepley ierr = DMPlexInsertBoundaryValues(plex, PETSC_TRUE, locX, time, faceGeometryFVM, NULL, NULL);CHKERRQ(ierr); 122a40652d4SToby Isaac /* TODO: locX_t */ 123c5d70e09SMatthew G. Knepley ierr = DMDestroy(&plex);CHKERRQ(ierr); 124c5d70e09SMatthew G. Knepley PetscFunctionReturn(0); 125c5d70e09SMatthew G. Knepley } 126c5d70e09SMatthew G. Knepley 12724cdb843SMatthew G. Knepley /*@ 12824cdb843SMatthew G. Knepley DMPlexTSComputeIFunctionFEM - Form the local residual F from the local input X using pointwise functions specified by the user 12924cdb843SMatthew G. Knepley 13024cdb843SMatthew G. Knepley Input Parameters: 13124cdb843SMatthew G. Knepley + dm - The mesh 13224cdb843SMatthew G. Knepley . t - The time 13308449791SMatthew G. Knepley . locX - Local solution 13408449791SMatthew G. Knepley . locX_t - Local solution time derivative, or NULL 13524cdb843SMatthew G. Knepley - user - The user context 13624cdb843SMatthew G. Knepley 13724cdb843SMatthew G. Knepley Output Parameter: 13808449791SMatthew G. Knepley . locF - Local output vector 13924cdb843SMatthew G. Knepley 14024cdb843SMatthew G. Knepley Level: developer 14124cdb843SMatthew G. Knepley 14224cdb843SMatthew G. Knepley .seealso: DMPlexComputeJacobianActionFEM() 14324cdb843SMatthew G. Knepley @*/ 14408449791SMatthew G. Knepley PetscErrorCode DMPlexTSComputeIFunctionFEM(DM dm, PetscReal time, Vec locX, Vec locX_t, Vec locF, void *user) 14524cdb843SMatthew G. Knepley { 1466da023fcSToby Isaac DM plex; 1474a3e9fdbSToby Isaac IS cellIS; 1484a3e9fdbSToby Isaac PetscInt depth; 14924cdb843SMatthew G. Knepley PetscErrorCode ierr; 15024cdb843SMatthew G. Knepley 15124cdb843SMatthew G. Knepley PetscFunctionBegin; 1526da023fcSToby Isaac ierr = DMTSConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr); 1534a3e9fdbSToby Isaac ierr = DMPlexGetDepth(plex, &depth);CHKERRQ(ierr); 154aeadca18SToby Isaac ierr = DMGetStratumIS(plex, "dim", depth, &cellIS);CHKERRQ(ierr); 1554a3e9fdbSToby Isaac if (!cellIS) { 1564a3e9fdbSToby Isaac ierr = DMGetStratumIS(plex, "depth", depth, &cellIS);CHKERRQ(ierr); 1574a3e9fdbSToby Isaac } 1584a3e9fdbSToby Isaac ierr = DMPlexComputeResidual_Internal(plex, cellIS, time, locX, locX_t, time, locF, user);CHKERRQ(ierr); 1594a3e9fdbSToby Isaac ierr = ISDestroy(&cellIS);CHKERRQ(ierr); 160b2338dbcSToby Isaac ierr = DMDestroy(&plex);CHKERRQ(ierr); 16124cdb843SMatthew G. Knepley PetscFunctionReturn(0); 16224cdb843SMatthew G. Knepley } 1637cdb2a12SMatthew G. Knepley 164756a1f44SMatthew G. Knepley /*@ 165756a1f44SMatthew G. Knepley DMPlexTSComputeIJacobianFEM - Form the local Jacobian J from the local input X using pointwise functions specified by the user 166756a1f44SMatthew G. Knepley 167756a1f44SMatthew G. Knepley Input Parameters: 168756a1f44SMatthew G. Knepley + dm - The mesh 169756a1f44SMatthew G. Knepley . t - The time 170756a1f44SMatthew G. Knepley . locX - Local solution 171756a1f44SMatthew G. Knepley . locX_t - Local solution time derivative, or NULL 172756a1f44SMatthew G. Knepley . X_tshift - The multiplicative parameter for dF/du_t 173756a1f44SMatthew G. Knepley - user - The user context 174756a1f44SMatthew G. Knepley 175756a1f44SMatthew G. Knepley Output Parameter: 176756a1f44SMatthew G. Knepley . locF - Local output vector 177756a1f44SMatthew G. Knepley 178756a1f44SMatthew G. Knepley Level: developer 179756a1f44SMatthew G. Knepley 180756a1f44SMatthew G. Knepley .seealso: DMPlexComputeJacobianActionFEM() 181756a1f44SMatthew G. Knepley @*/ 182756a1f44SMatthew G. Knepley PetscErrorCode DMPlexTSComputeIJacobianFEM(DM dm, PetscReal time, Vec locX, Vec locX_t, PetscReal X_tShift, Mat Jac, Mat JacP, void *user) 183756a1f44SMatthew G. Knepley { 184756a1f44SMatthew G. Knepley DM plex; 185f04eb4edSMatthew G. Knepley PetscDS prob; 186f04eb4edSMatthew G. Knepley PetscBool hasJac, hasPrec; 1879044fa66SMatthew G. Knepley IS cellIS; 1884a3e9fdbSToby Isaac PetscInt depth; 189756a1f44SMatthew G. Knepley PetscErrorCode ierr; 190756a1f44SMatthew G. Knepley 191756a1f44SMatthew G. Knepley PetscFunctionBegin; 192756a1f44SMatthew G. Knepley ierr = DMTSConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr); 193f04eb4edSMatthew G. Knepley ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 194f04eb4edSMatthew G. Knepley ierr = PetscDSHasJacobian(prob, &hasJac);CHKERRQ(ierr); 195f04eb4edSMatthew G. Knepley ierr = PetscDSHasJacobianPreconditioner(prob, &hasPrec);CHKERRQ(ierr); 196f04eb4edSMatthew G. Knepley if (hasJac && hasPrec) {ierr = MatZeroEntries(Jac);CHKERRQ(ierr);} 197f04eb4edSMatthew G. Knepley ierr = MatZeroEntries(JacP);CHKERRQ(ierr); 1984a3e9fdbSToby Isaac ierr = DMPlexGetDepth(plex,&depth);CHKERRQ(ierr); 199aeadca18SToby Isaac ierr = DMGetStratumIS(plex, "dim", depth, &cellIS);CHKERRQ(ierr); 2009044fa66SMatthew G. Knepley if (!cellIS) {ierr = DMGetStratumIS(plex, "depth", depth, &cellIS);CHKERRQ(ierr);} 2014a3e9fdbSToby Isaac ierr = DMPlexComputeJacobian_Internal(plex, cellIS, time, X_tShift, locX, locX_t, Jac, JacP, user);CHKERRQ(ierr); 2024a3e9fdbSToby Isaac ierr = ISDestroy(&cellIS);CHKERRQ(ierr); 203756a1f44SMatthew G. Knepley ierr = DMDestroy(&plex);CHKERRQ(ierr); 204756a1f44SMatthew G. Knepley PetscFunctionReturn(0); 205756a1f44SMatthew G. Knepley } 206756a1f44SMatthew G. Knepley 207bee9a294SMatthew G. Knepley /*@C 208*f2cacb80SMatthew G. Knepley DMTSCheckResidual - Check the residual of the exact solution 209*f2cacb80SMatthew G. Knepley 210*f2cacb80SMatthew G. Knepley Input Parameters: 211*f2cacb80SMatthew G. Knepley + ts - the TS object 212*f2cacb80SMatthew G. Knepley . dm - the DM 213*f2cacb80SMatthew G. Knepley . t - the time 214*f2cacb80SMatthew G. Knepley . u - a DM vector 215*f2cacb80SMatthew G. Knepley . u_t - a DM vector 216*f2cacb80SMatthew G. Knepley - tol - A tolerance for the check, or -1 to print the results instead 217*f2cacb80SMatthew G. Knepley 218*f2cacb80SMatthew G. Knepley Output Parameters: 219*f2cacb80SMatthew G. Knepley . residual - The residual norm of the exact solution, or NULL 220*f2cacb80SMatthew G. Knepley 221*f2cacb80SMatthew G. Knepley Level: developer 222*f2cacb80SMatthew G. Knepley 223*f2cacb80SMatthew G. Knepley .seealso: DNTSCheckFromOptions(), DMTSCheckJacobian(), DNSNESCheckFromOptions(), DMSNESCheckDiscretization(), DMSNESCheckJacobian() 224*f2cacb80SMatthew G. Knepley @*/ 225*f2cacb80SMatthew G. Knepley PetscErrorCode DMTSCheckResidual(TS ts, DM dm, PetscReal t, Vec u, Vec u_t, PetscReal tol, PetscReal *residual) 226*f2cacb80SMatthew G. Knepley { 227*f2cacb80SMatthew G. Knepley MPI_Comm comm; 228*f2cacb80SMatthew G. Knepley Vec r; 229*f2cacb80SMatthew G. Knepley PetscReal res; 230*f2cacb80SMatthew G. Knepley PetscErrorCode ierr; 231*f2cacb80SMatthew G. Knepley 232*f2cacb80SMatthew G. Knepley PetscFunctionBegin; 233*f2cacb80SMatthew G. Knepley PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 234*f2cacb80SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 2); 235*f2cacb80SMatthew G. Knepley PetscValidHeaderSpecific(u, VEC_CLASSID, 3); 236*f2cacb80SMatthew G. Knepley if (residual) PetscValidRealPointer(residual, 5); 237*f2cacb80SMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject) ts, &comm);CHKERRQ(ierr); 238*f2cacb80SMatthew G. Knepley ierr = DMComputeExactSolution(dm, t, u, u_t);CHKERRQ(ierr); 239*f2cacb80SMatthew G. Knepley ierr = VecDuplicate(u, &r);CHKERRQ(ierr); 240*f2cacb80SMatthew G. Knepley ierr = TSComputeIFunction(ts, t, u, u_t, r, PETSC_FALSE);CHKERRQ(ierr); 241*f2cacb80SMatthew G. Knepley ierr = VecNorm(r, NORM_2, &res);CHKERRQ(ierr); 242*f2cacb80SMatthew G. Knepley if (tol >= 0.0) { 243*f2cacb80SMatthew G. Knepley if (res > tol) SETERRQ2(comm, PETSC_ERR_ARG_WRONG, "L_2 Residual %g exceeds tolerance %g", (double) res, (double) tol); 244*f2cacb80SMatthew G. Knepley } else if (residual) { 245*f2cacb80SMatthew G. Knepley *residual = res; 246*f2cacb80SMatthew G. Knepley } else { 247*f2cacb80SMatthew G. Knepley ierr = PetscPrintf(comm, "L_2 Residual: %g\n", (double)res);CHKERRQ(ierr); 248*f2cacb80SMatthew G. Knepley ierr = VecChop(r, 1.0e-10);CHKERRQ(ierr); 249*f2cacb80SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) r, "Initial Residual");CHKERRQ(ierr); 250*f2cacb80SMatthew G. Knepley ierr = PetscObjectSetOptionsPrefix((PetscObject)r,"res_");CHKERRQ(ierr); 251*f2cacb80SMatthew G. Knepley ierr = VecViewFromOptions(r, NULL, "-vec_view");CHKERRQ(ierr); 252*f2cacb80SMatthew G. Knepley } 253*f2cacb80SMatthew G. Knepley ierr = VecDestroy(&r);CHKERRQ(ierr); 254*f2cacb80SMatthew G. Knepley PetscFunctionReturn(0); 255*f2cacb80SMatthew G. Knepley } 256*f2cacb80SMatthew G. Knepley 257*f2cacb80SMatthew G. Knepley /*@C 258*f2cacb80SMatthew G. Knepley DMTSCheckJacobian - Check the Jacobian of the exact solution against the residual using the Taylor Test 259*f2cacb80SMatthew G. Knepley 260*f2cacb80SMatthew G. Knepley Input Parameters: 261*f2cacb80SMatthew G. Knepley + ts - the TS object 262*f2cacb80SMatthew G. Knepley . dm - the DM 263*f2cacb80SMatthew G. Knepley . t - the time 264*f2cacb80SMatthew G. Knepley . u - a DM vector 265*f2cacb80SMatthew G. Knepley . u_t - a DM vector 266*f2cacb80SMatthew G. Knepley - tol - A tolerance for the check, or -1 to print the results instead 267*f2cacb80SMatthew G. Knepley 268*f2cacb80SMatthew G. Knepley Output Parameters: 269*f2cacb80SMatthew G. Knepley + isLinear - Flag indicaing that the function looks linear, or NULL 270*f2cacb80SMatthew G. Knepley - convRate - The rate of convergence of the linear model, or NULL 271*f2cacb80SMatthew G. Knepley 272*f2cacb80SMatthew G. Knepley Level: developer 273*f2cacb80SMatthew G. Knepley 274*f2cacb80SMatthew G. Knepley .seealso: DNTSCheckFromOptions(), DMTSCheckResidual(), DNSNESCheckFromOptions(), DMSNESCheckDiscretization(), DMSNESCheckResidual() 275*f2cacb80SMatthew G. Knepley @*/ 276*f2cacb80SMatthew G. Knepley PetscErrorCode DMTSCheckJacobian(TS ts, DM dm, PetscReal t, Vec u, Vec u_t, PetscReal tol, PetscBool *isLinear, PetscReal *convRate) 277*f2cacb80SMatthew G. Knepley { 278*f2cacb80SMatthew G. Knepley MPI_Comm comm; 279*f2cacb80SMatthew G. Knepley PetscDS ds; 280*f2cacb80SMatthew G. Knepley Mat J, M; 281*f2cacb80SMatthew G. Knepley MatNullSpace nullspace; 282*f2cacb80SMatthew G. Knepley PetscReal dt, shift, slope, intercept; 283*f2cacb80SMatthew G. Knepley PetscBool hasJac, hasPrec, isLin = PETSC_FALSE; 284*f2cacb80SMatthew G. Knepley PetscErrorCode ierr; 285*f2cacb80SMatthew G. Knepley 286*f2cacb80SMatthew G. Knepley PetscFunctionBegin; 287*f2cacb80SMatthew G. Knepley PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 288*f2cacb80SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 2); 289*f2cacb80SMatthew G. Knepley PetscValidHeaderSpecific(u, VEC_CLASSID, 3); 290*f2cacb80SMatthew G. Knepley if (isLinear) PetscValidBoolPointer(isLinear, 5); 291*f2cacb80SMatthew G. Knepley if (convRate) PetscValidRealPointer(convRate, 5); 292*f2cacb80SMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject) ts, &comm);CHKERRQ(ierr); 293*f2cacb80SMatthew G. Knepley ierr = DMComputeExactSolution(dm, t, u, u_t);CHKERRQ(ierr); 294*f2cacb80SMatthew G. Knepley /* Create and view matrices */ 295*f2cacb80SMatthew G. Knepley ierr = TSGetTimeStep(ts, &dt);CHKERRQ(ierr); 296*f2cacb80SMatthew G. Knepley shift = 1.0/dt; 297*f2cacb80SMatthew G. Knepley ierr = DMCreateMatrix(dm, &J);CHKERRQ(ierr); 298*f2cacb80SMatthew G. Knepley ierr = DMGetDS(dm, &ds);CHKERRQ(ierr); 299*f2cacb80SMatthew G. Knepley ierr = PetscDSHasJacobian(ds, &hasJac);CHKERRQ(ierr); 300*f2cacb80SMatthew G. Knepley ierr = PetscDSHasJacobianPreconditioner(ds, &hasPrec);CHKERRQ(ierr); 301*f2cacb80SMatthew G. Knepley if (hasJac && hasPrec) { 302*f2cacb80SMatthew G. Knepley ierr = DMCreateMatrix(dm, &M);CHKERRQ(ierr); 303*f2cacb80SMatthew G. Knepley ierr = TSComputeIJacobian(ts, t, u, u_t, shift, J, M, PETSC_FALSE);CHKERRQ(ierr); 304*f2cacb80SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) M, "Preconditioning Matrix");CHKERRQ(ierr); 305*f2cacb80SMatthew G. Knepley ierr = PetscObjectSetOptionsPrefix((PetscObject) M, "jacpre_");CHKERRQ(ierr); 306*f2cacb80SMatthew G. Knepley ierr = MatViewFromOptions(M, NULL, "-mat_view");CHKERRQ(ierr); 307*f2cacb80SMatthew G. Knepley ierr = MatDestroy(&M);CHKERRQ(ierr); 308*f2cacb80SMatthew G. Knepley } else { 309*f2cacb80SMatthew G. Knepley ierr = TSComputeIJacobian(ts, t, u, u_t, shift, J, J, PETSC_FALSE);CHKERRQ(ierr); 310*f2cacb80SMatthew G. Knepley } 311*f2cacb80SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) J, "Jacobian");CHKERRQ(ierr); 312*f2cacb80SMatthew G. Knepley ierr = PetscObjectSetOptionsPrefix((PetscObject) J, "jac_");CHKERRQ(ierr); 313*f2cacb80SMatthew G. Knepley ierr = MatViewFromOptions(J, NULL, "-mat_view");CHKERRQ(ierr); 314*f2cacb80SMatthew G. Knepley /* Check nullspace */ 315*f2cacb80SMatthew G. Knepley ierr = MatGetNullSpace(J, &nullspace);CHKERRQ(ierr); 316*f2cacb80SMatthew G. Knepley if (nullspace) { 317*f2cacb80SMatthew G. Knepley PetscBool isNull; 318*f2cacb80SMatthew G. Knepley ierr = MatNullSpaceTest(nullspace, J, &isNull);CHKERRQ(ierr); 319*f2cacb80SMatthew G. Knepley if (!isNull) SETERRQ(comm, PETSC_ERR_PLIB, "The null space calculated for the system operator is invalid."); 320*f2cacb80SMatthew G. Knepley } 321*f2cacb80SMatthew G. Knepley ierr = MatNullSpaceDestroy(&nullspace);CHKERRQ(ierr); 322*f2cacb80SMatthew G. Knepley /* Taylor test */ 323*f2cacb80SMatthew G. Knepley { 324*f2cacb80SMatthew G. Knepley PetscRandom rand; 325*f2cacb80SMatthew G. Knepley Vec du, uhat, uhat_t, r, rhat, df; 326*f2cacb80SMatthew G. Knepley PetscReal h; 327*f2cacb80SMatthew G. Knepley PetscReal *es, *hs, *errors; 328*f2cacb80SMatthew G. Knepley PetscReal hMax = 1.0, hMin = 1e-6, hMult = 0.1; 329*f2cacb80SMatthew G. Knepley PetscInt Nv, v; 330*f2cacb80SMatthew G. Knepley 331*f2cacb80SMatthew G. Knepley /* Choose a perturbation direction */ 332*f2cacb80SMatthew G. Knepley ierr = PetscRandomCreate(comm, &rand);CHKERRQ(ierr); 333*f2cacb80SMatthew G. Knepley ierr = VecDuplicate(u, &du);CHKERRQ(ierr); 334*f2cacb80SMatthew G. Knepley ierr = VecSetRandom(du, rand); CHKERRQ(ierr); 335*f2cacb80SMatthew G. Knepley ierr = PetscRandomDestroy(&rand);CHKERRQ(ierr); 336*f2cacb80SMatthew G. Knepley ierr = VecDuplicate(u, &df);CHKERRQ(ierr); 337*f2cacb80SMatthew G. Knepley ierr = MatMult(J, du, df);CHKERRQ(ierr); 338*f2cacb80SMatthew G. Knepley /* Evaluate residual at u, F(u), save in vector r */ 339*f2cacb80SMatthew G. Knepley ierr = VecDuplicate(u, &r);CHKERRQ(ierr); 340*f2cacb80SMatthew G. Knepley ierr = TSComputeIFunction(ts, t, u, u_t, r, PETSC_FALSE);CHKERRQ(ierr); 341*f2cacb80SMatthew G. Knepley /* Look at the convergence of our Taylor approximation as we approach u */ 342*f2cacb80SMatthew G. Knepley for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv); 343*f2cacb80SMatthew G. Knepley ierr = PetscCalloc3(Nv, &es, Nv, &hs, Nv, &errors);CHKERRQ(ierr); 344*f2cacb80SMatthew G. Knepley ierr = VecDuplicate(u, &uhat);CHKERRQ(ierr); 345*f2cacb80SMatthew G. Knepley ierr = VecDuplicate(u, &uhat_t);CHKERRQ(ierr); 346*f2cacb80SMatthew G. Knepley ierr = VecDuplicate(u, &rhat);CHKERRQ(ierr); 347*f2cacb80SMatthew G. Knepley for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv) { 348*f2cacb80SMatthew G. Knepley ierr = VecWAXPY(uhat, h, du, u);CHKERRQ(ierr); 349*f2cacb80SMatthew G. Knepley ierr = VecWAXPY(uhat_t, h*shift, du, u_t);CHKERRQ(ierr); 350*f2cacb80SMatthew G. Knepley /* F(\hat u, \hat u_t) \approx F(u, u_t) + J(u, u_t) (uhat - u) + J_t(u, u_t) (uhat_t - u_t) = F(u) + h * J(u) du + h * shift * J_t(u) du = F(u) + h F' du */ 351*f2cacb80SMatthew G. Knepley ierr = TSComputeIFunction(ts, t, uhat, uhat_t, rhat, PETSC_FALSE);CHKERRQ(ierr); 352*f2cacb80SMatthew G. Knepley ierr = VecAXPBYPCZ(rhat, -1.0, -h, 1.0, r, df);CHKERRQ(ierr); 353*f2cacb80SMatthew G. Knepley ierr = VecNorm(rhat, NORM_2, &errors[Nv]);CHKERRQ(ierr); 354*f2cacb80SMatthew G. Knepley 355*f2cacb80SMatthew G. Knepley es[Nv] = PetscLog10Real(errors[Nv]); 356*f2cacb80SMatthew G. Knepley hs[Nv] = PetscLog10Real(h); 357*f2cacb80SMatthew G. Knepley } 358*f2cacb80SMatthew G. Knepley ierr = VecDestroy(&uhat);CHKERRQ(ierr); 359*f2cacb80SMatthew G. Knepley ierr = VecDestroy(&uhat_t);CHKERRQ(ierr); 360*f2cacb80SMatthew G. Knepley ierr = VecDestroy(&rhat);CHKERRQ(ierr); 361*f2cacb80SMatthew G. Knepley ierr = VecDestroy(&df);CHKERRQ(ierr); 362*f2cacb80SMatthew G. Knepley ierr = VecDestroy(&r);CHKERRQ(ierr); 363*f2cacb80SMatthew G. Knepley ierr = VecDestroy(&du);CHKERRQ(ierr); 364*f2cacb80SMatthew G. Knepley for (v = 0; v < Nv; ++v) { 365*f2cacb80SMatthew G. Knepley if ((tol >= 0) && (errors[v] > tol)) break; 366*f2cacb80SMatthew G. Knepley else if (errors[v] > PETSC_SMALL) break; 367*f2cacb80SMatthew G. Knepley } 368*f2cacb80SMatthew G. Knepley if (v == Nv) isLin = PETSC_TRUE; 369*f2cacb80SMatthew G. Knepley ierr = PetscLinearRegression(Nv, hs, es, &slope, &intercept);CHKERRQ(ierr); 370*f2cacb80SMatthew G. Knepley ierr = PetscFree3(es, hs, errors);CHKERRQ(ierr); 371*f2cacb80SMatthew G. Knepley /* Slope should be about 2 */ 372*f2cacb80SMatthew G. Knepley if (tol >= 0) { 373*f2cacb80SMatthew G. Knepley if (!isLin && PetscAbsReal(2 - slope) > tol) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Taylor approximation convergence rate should be 2, not %0.2f", (double) slope); 374*f2cacb80SMatthew G. Knepley } else if (isLinear || convRate) { 375*f2cacb80SMatthew G. Knepley if (isLinear) *isLinear = isLin; 376*f2cacb80SMatthew G. Knepley if (convRate) *convRate = slope; 377*f2cacb80SMatthew G. Knepley } else { 378*f2cacb80SMatthew G. Knepley if (!isLin) {ierr = PetscPrintf(comm, "Taylor approximation converging at order %3.2f\n", (double) slope);CHKERRQ(ierr);} 379*f2cacb80SMatthew G. Knepley else {ierr = PetscPrintf(comm, "Function appears to be linear\n");CHKERRQ(ierr);} 380*f2cacb80SMatthew G. Knepley } 381*f2cacb80SMatthew G. Knepley } 382*f2cacb80SMatthew G. Knepley ierr = MatDestroy(&J);CHKERRQ(ierr); 383*f2cacb80SMatthew G. Knepley PetscFunctionReturn(0); 384*f2cacb80SMatthew G. Knepley } 385*f2cacb80SMatthew G. Knepley 386*f2cacb80SMatthew G. Knepley /*@C 387bee9a294SMatthew G. Knepley DMTSCheckFromOptions - Check the residual and Jacobian functions using the exact solution by outputting some diagnostic information 388bee9a294SMatthew G. Knepley 389bee9a294SMatthew G. Knepley Input Parameters: 390bee9a294SMatthew G. Knepley + ts - the TS object 3917f96f943SMatthew G. Knepley - u - representative TS vector 3927f96f943SMatthew G. Knepley 3937f96f943SMatthew G. Knepley Note: The user must call PetscDSSetExactSolution() beforehand 394bee9a294SMatthew G. Knepley 395bee9a294SMatthew G. Knepley Level: developer 396bee9a294SMatthew G. Knepley @*/ 3977f96f943SMatthew G. Knepley PetscErrorCode DMTSCheckFromOptions(TS ts, Vec u) 3987cdb2a12SMatthew G. Knepley { 3997cdb2a12SMatthew G. Knepley DM dm; 4007cdb2a12SMatthew G. Knepley SNES snes; 401*f2cacb80SMatthew G. Knepley Vec sol, u_t; 402*f2cacb80SMatthew G. Knepley PetscReal t; 4037cdb2a12SMatthew G. Knepley PetscBool check; 4047cdb2a12SMatthew G. Knepley PetscErrorCode ierr; 4057cdb2a12SMatthew G. Knepley 4067cdb2a12SMatthew G. Knepley PetscFunctionBegin; 407c5929fdfSBarry Smith ierr = PetscOptionsHasName(((PetscObject)ts)->options,((PetscObject)ts)->prefix, "-dmts_check", &check);CHKERRQ(ierr); 4087cdb2a12SMatthew G. Knepley if (!check) PetscFunctionReturn(0); 40972fd7fc8SMatthew G. Knepley ierr = VecDuplicate(u, &sol);CHKERRQ(ierr); 4106fb91e28SMatthew G. Knepley ierr = VecCopy(u, sol);CHKERRQ(ierr); 4116fb91e28SMatthew G. Knepley ierr = TSSetSolution(ts, u);CHKERRQ(ierr); 4127cdb2a12SMatthew G. Knepley ierr = TSGetDM(ts, &dm);CHKERRQ(ierr); 4137cdb2a12SMatthew G. Knepley ierr = TSSetUp(ts);CHKERRQ(ierr); 4147cdb2a12SMatthew G. Knepley ierr = TSGetSNES(ts, &snes);CHKERRQ(ierr); 4156fb91e28SMatthew G. Knepley ierr = SNESSetSolution(snes, u);CHKERRQ(ierr); 416*f2cacb80SMatthew G. Knepley 417*f2cacb80SMatthew G. Knepley ierr = TSGetTime(ts, &t);CHKERRQ(ierr); 418*f2cacb80SMatthew G. Knepley ierr = DMSNESCheckDiscretization(snes, dm, t, sol, -1.0, NULL);CHKERRQ(ierr); 419*f2cacb80SMatthew G. Knepley ierr = DMGetGlobalVector(dm, &u_t);CHKERRQ(ierr); 420*f2cacb80SMatthew G. Knepley ierr = DMTSCheckResidual(ts, dm, t, sol, u_t, -1.0, NULL);CHKERRQ(ierr); 421*f2cacb80SMatthew G. Knepley ierr = DMTSCheckJacobian(ts, dm, t, sol, u_t, -1.0, NULL, NULL);CHKERRQ(ierr); 422*f2cacb80SMatthew G. Knepley ierr = DMRestoreGlobalVector(dm, &u_t);CHKERRQ(ierr); 423*f2cacb80SMatthew G. Knepley 42472fd7fc8SMatthew G. Knepley ierr = VecDestroy(&sol);CHKERRQ(ierr); 4257cdb2a12SMatthew G. Knepley PetscFunctionReturn(0); 4267cdb2a12SMatthew G. Knepley } 427