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; 63*6528b96dSMatthew G. Knepley PetscHashFormKey key = {NULL, 0, 0}; 64254c1ad2SMatthew G. Knepley PetscErrorCode ierr; 65254c1ad2SMatthew G. Knepley 66254c1ad2SMatthew G. Knepley PetscFunctionBegin; 676da023fcSToby Isaac ierr = DMTSConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr); 684a3e9fdbSToby Isaac ierr = DMPlexGetDepth(plex, &depth);CHKERRQ(ierr); 69aeadca18SToby Isaac ierr = DMGetStratumIS(plex, "dim", depth, &cellIS);CHKERRQ(ierr); 704a3e9fdbSToby Isaac if (!cellIS) { 714a3e9fdbSToby Isaac ierr = DMGetStratumIS(plex, "depth", depth, &cellIS);CHKERRQ(ierr); 724a3e9fdbSToby Isaac } 736da023fcSToby Isaac ierr = DMGetLocalVector(plex, &locF);CHKERRQ(ierr); 743b16df42SMatthew G. Knepley ierr = VecZeroEntries(locF);CHKERRQ(ierr); 75*6528b96dSMatthew G. Knepley ierr = DMPlexComputeResidual_Internal(plex, key, cellIS, time, locX, NULL, time, locF, user);CHKERRQ(ierr); 769bda831aSToby Isaac ierr = DMLocalToGlobalBegin(plex, locF, ADD_VALUES, F);CHKERRQ(ierr); 779bda831aSToby Isaac ierr = DMLocalToGlobalEnd(plex, locF, ADD_VALUES, F);CHKERRQ(ierr); 786da023fcSToby Isaac ierr = DMRestoreLocalVector(plex, &locF);CHKERRQ(ierr); 794a3e9fdbSToby Isaac ierr = ISDestroy(&cellIS);CHKERRQ(ierr); 80b2338dbcSToby Isaac ierr = DMDestroy(&plex);CHKERRQ(ierr); 81254c1ad2SMatthew G. Knepley PetscFunctionReturn(0); 82254c1ad2SMatthew G. Knepley } 83254c1ad2SMatthew G. Knepley 84c5d70e09SMatthew G. Knepley /*@ 85ef68eab9SMatthew 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 86c5d70e09SMatthew G. Knepley 87c5d70e09SMatthew G. Knepley Input Parameters: 88c5d70e09SMatthew G. Knepley + dm - The mesh 89a40652d4SToby Isaac . t - The time 90a40652d4SToby Isaac . locX - Local solution 91a40652d4SToby Isaac . locX_t - Local solution time derivative, or NULL 92c5d70e09SMatthew G. Knepley - user - The user context 93c5d70e09SMatthew G. Knepley 94c5d70e09SMatthew G. Knepley Level: developer 95c5d70e09SMatthew G. Knepley 96a40652d4SToby Isaac .seealso: DMPlexComputeJacobianActionFEM() 97c5d70e09SMatthew G. Knepley @*/ 98ef68eab9SMatthew G. Knepley PetscErrorCode DMPlexTSComputeBoundary(DM dm, PetscReal time, Vec locX, Vec locX_t, void *user) 99c5d70e09SMatthew G. Knepley { 100c5d70e09SMatthew G. Knepley DM plex; 101ef68eab9SMatthew G. Knepley Vec faceGeometryFVM = NULL; 102ef68eab9SMatthew G. Knepley PetscInt Nf, f; 103c5d70e09SMatthew G. Knepley PetscErrorCode ierr; 104c5d70e09SMatthew G. Knepley 105c5d70e09SMatthew G. Knepley PetscFunctionBegin; 106c5d70e09SMatthew G. Knepley ierr = DMTSConvertPlex(dm, &plex, PETSC_TRUE);CHKERRQ(ierr); 107ef68eab9SMatthew G. Knepley ierr = DMGetNumFields(plex, &Nf);CHKERRQ(ierr); 1089586001cSMatthew G. Knepley if (!locX_t) { 1099586001cSMatthew G. Knepley /* This is the RHS part */ 110ef68eab9SMatthew G. Knepley for (f = 0; f < Nf; f++) { 111ef68eab9SMatthew G. Knepley PetscObject obj; 112ef68eab9SMatthew G. Knepley PetscClassId id; 113ef68eab9SMatthew G. Knepley 11444a7f3ddSMatthew G. Knepley ierr = DMGetField(plex, f, NULL, &obj);CHKERRQ(ierr); 115ef68eab9SMatthew G. Knepley ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 116ef68eab9SMatthew G. Knepley if (id == PETSCFV_CLASSID) { 1173e9753d6SMatthew G. Knepley ierr = DMPlexGetGeometryFVM(plex, &faceGeometryFVM, NULL, NULL);CHKERRQ(ierr); 118ef68eab9SMatthew G. Knepley break; 119ef68eab9SMatthew G. Knepley } 120ef68eab9SMatthew G. Knepley } 1219586001cSMatthew G. Knepley } 1229586001cSMatthew G. Knepley ierr = DMPlexInsertBoundaryValues(plex, PETSC_TRUE, locX, time, faceGeometryFVM, NULL, NULL);CHKERRQ(ierr); 12356cf3b9cSMatthew G. Knepley ierr = DMPlexInsertTimeDerivativeBoundaryValues(plex, PETSC_TRUE, locX_t, time, faceGeometryFVM, NULL, NULL);CHKERRQ(ierr); 124c5d70e09SMatthew G. Knepley ierr = DMDestroy(&plex);CHKERRQ(ierr); 125c5d70e09SMatthew G. Knepley PetscFunctionReturn(0); 126c5d70e09SMatthew G. Knepley } 127c5d70e09SMatthew G. Knepley 12824cdb843SMatthew G. Knepley /*@ 12924cdb843SMatthew G. Knepley DMPlexTSComputeIFunctionFEM - Form the local residual F from the local input X using pointwise functions specified by the user 13024cdb843SMatthew G. Knepley 13124cdb843SMatthew G. Knepley Input Parameters: 13224cdb843SMatthew G. Knepley + dm - The mesh 13324cdb843SMatthew G. Knepley . t - The time 13408449791SMatthew G. Knepley . locX - Local solution 13508449791SMatthew G. Knepley . locX_t - Local solution time derivative, or NULL 13624cdb843SMatthew G. Knepley - user - The user context 13724cdb843SMatthew G. Knepley 13824cdb843SMatthew G. Knepley Output Parameter: 13908449791SMatthew G. Knepley . locF - Local output vector 14024cdb843SMatthew G. Knepley 14124cdb843SMatthew G. Knepley Level: developer 14224cdb843SMatthew G. Knepley 14324cdb843SMatthew G. Knepley .seealso: DMPlexComputeJacobianActionFEM() 14424cdb843SMatthew G. Knepley @*/ 14508449791SMatthew G. Knepley PetscErrorCode DMPlexTSComputeIFunctionFEM(DM dm, PetscReal time, Vec locX, Vec locX_t, Vec locF, void *user) 14624cdb843SMatthew G. Knepley { 1476da023fcSToby Isaac DM plex; 148*6528b96dSMatthew G. Knepley IS allcellIS; 149*6528b96dSMatthew G. Knepley PetscInt Nds, s; 15024cdb843SMatthew G. Knepley PetscErrorCode ierr; 15124cdb843SMatthew G. Knepley 15224cdb843SMatthew G. Knepley PetscFunctionBegin; 1536da023fcSToby Isaac ierr = DMTSConvertPlex(dm, &plex, PETSC_TRUE);CHKERRQ(ierr); 154*6528b96dSMatthew G. Knepley ierr = DMPlexGetAllCells_Internal(plex, &allcellIS);CHKERRQ(ierr); 155*6528b96dSMatthew G. Knepley ierr = DMGetNumDS(dm, &Nds);CHKERRQ(ierr); 156*6528b96dSMatthew G. Knepley for (s = 0; s < Nds; ++s) { 157*6528b96dSMatthew G. Knepley PetscDS ds; 158*6528b96dSMatthew G. Knepley IS cellIS; 159*6528b96dSMatthew G. Knepley PetscHashFormKey key; 160*6528b96dSMatthew G. Knepley 161*6528b96dSMatthew G. Knepley ierr = DMGetRegionNumDS(dm, s, &key.label, NULL, &ds);CHKERRQ(ierr); 162*6528b96dSMatthew G. Knepley key.value = 0; 163*6528b96dSMatthew G. Knepley key.field = 0; 164*6528b96dSMatthew G. Knepley if (!key.label) { 165*6528b96dSMatthew G. Knepley ierr = PetscObjectReference((PetscObject) allcellIS);CHKERRQ(ierr); 166*6528b96dSMatthew G. Knepley cellIS = allcellIS; 167*6528b96dSMatthew G. Knepley } else { 168*6528b96dSMatthew G. Knepley IS pointIS; 169*6528b96dSMatthew G. Knepley 170*6528b96dSMatthew G. Knepley key.value = 1; 171*6528b96dSMatthew G. Knepley ierr = DMLabelGetStratumIS(key.label, key.value, &pointIS);CHKERRQ(ierr); 172*6528b96dSMatthew G. Knepley ierr = ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS);CHKERRQ(ierr); 173*6528b96dSMatthew G. Knepley ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 1744a3e9fdbSToby Isaac } 175*6528b96dSMatthew G. Knepley ierr = DMPlexComputeResidual_Internal(plex, key, cellIS, time, locX, locX_t, time, locF, user);CHKERRQ(ierr); 1764a3e9fdbSToby Isaac ierr = ISDestroy(&cellIS);CHKERRQ(ierr); 177*6528b96dSMatthew G. Knepley } 178*6528b96dSMatthew G. Knepley ierr = ISDestroy(&allcellIS);CHKERRQ(ierr); 179b2338dbcSToby Isaac ierr = DMDestroy(&plex);CHKERRQ(ierr); 18024cdb843SMatthew G. Knepley PetscFunctionReturn(0); 18124cdb843SMatthew G. Knepley } 1827cdb2a12SMatthew G. Knepley 183756a1f44SMatthew G. Knepley /*@ 184756a1f44SMatthew G. Knepley DMPlexTSComputeIJacobianFEM - Form the local Jacobian J from the local input X using pointwise functions specified by the user 185756a1f44SMatthew G. Knepley 186756a1f44SMatthew G. Knepley Input Parameters: 187756a1f44SMatthew G. Knepley + dm - The mesh 188756a1f44SMatthew G. Knepley . t - The time 189756a1f44SMatthew G. Knepley . locX - Local solution 190756a1f44SMatthew G. Knepley . locX_t - Local solution time derivative, or NULL 191756a1f44SMatthew G. Knepley . X_tshift - The multiplicative parameter for dF/du_t 192756a1f44SMatthew G. Knepley - user - The user context 193756a1f44SMatthew G. Knepley 194756a1f44SMatthew G. Knepley Output Parameter: 195756a1f44SMatthew G. Knepley . locF - Local output vector 196756a1f44SMatthew G. Knepley 197756a1f44SMatthew G. Knepley Level: developer 198756a1f44SMatthew G. Knepley 199756a1f44SMatthew G. Knepley .seealso: DMPlexComputeJacobianActionFEM() 200756a1f44SMatthew G. Knepley @*/ 201756a1f44SMatthew G. Knepley PetscErrorCode DMPlexTSComputeIJacobianFEM(DM dm, PetscReal time, Vec locX, Vec locX_t, PetscReal X_tShift, Mat Jac, Mat JacP, void *user) 202756a1f44SMatthew G. Knepley { 203756a1f44SMatthew G. Knepley DM plex; 204*6528b96dSMatthew G. Knepley IS allcellIS; 205f04eb4edSMatthew G. Knepley PetscBool hasJac, hasPrec; 206*6528b96dSMatthew G. Knepley PetscInt Nds, s; 207756a1f44SMatthew G. Knepley PetscErrorCode ierr; 208756a1f44SMatthew G. Knepley 209756a1f44SMatthew G. Knepley PetscFunctionBegin; 210756a1f44SMatthew G. Knepley ierr = DMTSConvertPlex(dm, &plex, PETSC_TRUE);CHKERRQ(ierr); 211*6528b96dSMatthew G. Knepley ierr = DMPlexGetAllCells_Internal(plex, &allcellIS);CHKERRQ(ierr); 212*6528b96dSMatthew G. Knepley ierr = DMGetNumDS(dm, &Nds);CHKERRQ(ierr); 213*6528b96dSMatthew G. Knepley for (s = 0; s < Nds; ++s) { 214*6528b96dSMatthew G. Knepley PetscDS ds; 215*6528b96dSMatthew G. Knepley IS cellIS; 216*6528b96dSMatthew G. Knepley PetscHashFormKey key; 217*6528b96dSMatthew G. Knepley 218*6528b96dSMatthew G. Knepley ierr = DMGetRegionNumDS(dm, s, &key.label, NULL, &ds);CHKERRQ(ierr); 219*6528b96dSMatthew G. Knepley key.value = 0; 220*6528b96dSMatthew G. Knepley key.field = 0; 221*6528b96dSMatthew G. Knepley if (!key.label) { 222*6528b96dSMatthew G. Knepley ierr = PetscObjectReference((PetscObject) allcellIS);CHKERRQ(ierr); 223*6528b96dSMatthew G. Knepley cellIS = allcellIS; 224*6528b96dSMatthew G. Knepley } else { 225*6528b96dSMatthew G. Knepley IS pointIS; 226*6528b96dSMatthew G. Knepley 227*6528b96dSMatthew G. Knepley key.value = 1; 228*6528b96dSMatthew G. Knepley ierr = DMLabelGetStratumIS(key.label, key.value, &pointIS);CHKERRQ(ierr); 229*6528b96dSMatthew G. Knepley ierr = ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS);CHKERRQ(ierr); 230*6528b96dSMatthew G. Knepley ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 231*6528b96dSMatthew G. Knepley } 232*6528b96dSMatthew G. Knepley if (!s) { 233*6528b96dSMatthew G. Knepley ierr = PetscDSHasJacobian(ds, &hasJac);CHKERRQ(ierr); 234*6528b96dSMatthew G. Knepley ierr = PetscDSHasJacobianPreconditioner(ds, &hasPrec);CHKERRQ(ierr); 235f04eb4edSMatthew G. Knepley if (hasJac && hasPrec) {ierr = MatZeroEntries(Jac);CHKERRQ(ierr);} 236f04eb4edSMatthew G. Knepley ierr = MatZeroEntries(JacP);CHKERRQ(ierr); 237*6528b96dSMatthew G. Knepley } 238*6528b96dSMatthew G. Knepley ierr = DMPlexComputeJacobian_Internal(plex, key, cellIS, time, X_tShift, locX, locX_t, Jac, JacP, user);CHKERRQ(ierr); 2394a3e9fdbSToby Isaac ierr = ISDestroy(&cellIS);CHKERRQ(ierr); 240*6528b96dSMatthew G. Knepley } 241*6528b96dSMatthew G. Knepley ierr = ISDestroy(&allcellIS);CHKERRQ(ierr); 242756a1f44SMatthew G. Knepley ierr = DMDestroy(&plex);CHKERRQ(ierr); 243756a1f44SMatthew G. Knepley PetscFunctionReturn(0); 244756a1f44SMatthew G. Knepley } 245756a1f44SMatthew G. Knepley 246bee9a294SMatthew G. Knepley /*@C 247f2cacb80SMatthew G. Knepley DMTSCheckResidual - Check the residual of the exact solution 248f2cacb80SMatthew G. Knepley 249f2cacb80SMatthew G. Knepley Input Parameters: 250f2cacb80SMatthew G. Knepley + ts - the TS object 251f2cacb80SMatthew G. Knepley . dm - the DM 252f2cacb80SMatthew G. Knepley . t - the time 253f2cacb80SMatthew G. Knepley . u - a DM vector 254f2cacb80SMatthew G. Knepley . u_t - a DM vector 255f2cacb80SMatthew G. Knepley - tol - A tolerance for the check, or -1 to print the results instead 256f2cacb80SMatthew G. Knepley 257f2cacb80SMatthew G. Knepley Output Parameters: 258f2cacb80SMatthew G. Knepley . residual - The residual norm of the exact solution, or NULL 259f2cacb80SMatthew G. Knepley 260f2cacb80SMatthew G. Knepley Level: developer 261f2cacb80SMatthew G. Knepley 262f2cacb80SMatthew G. Knepley .seealso: DNTSCheckFromOptions(), DMTSCheckJacobian(), DNSNESCheckFromOptions(), DMSNESCheckDiscretization(), DMSNESCheckJacobian() 263f2cacb80SMatthew G. Knepley @*/ 264f2cacb80SMatthew G. Knepley PetscErrorCode DMTSCheckResidual(TS ts, DM dm, PetscReal t, Vec u, Vec u_t, PetscReal tol, PetscReal *residual) 265f2cacb80SMatthew G. Knepley { 266f2cacb80SMatthew G. Knepley MPI_Comm comm; 267f2cacb80SMatthew G. Knepley Vec r; 268f2cacb80SMatthew G. Knepley PetscReal res; 269f2cacb80SMatthew G. Knepley PetscErrorCode ierr; 270f2cacb80SMatthew G. Knepley 271f2cacb80SMatthew G. Knepley PetscFunctionBegin; 272f2cacb80SMatthew G. Knepley PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 273f2cacb80SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 2); 274f2cacb80SMatthew G. Knepley PetscValidHeaderSpecific(u, VEC_CLASSID, 3); 275f2cacb80SMatthew G. Knepley if (residual) PetscValidRealPointer(residual, 5); 276f2cacb80SMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject) ts, &comm);CHKERRQ(ierr); 277f2cacb80SMatthew G. Knepley ierr = DMComputeExactSolution(dm, t, u, u_t);CHKERRQ(ierr); 278f2cacb80SMatthew G. Knepley ierr = VecDuplicate(u, &r);CHKERRQ(ierr); 279f2cacb80SMatthew G. Knepley ierr = TSComputeIFunction(ts, t, u, u_t, r, PETSC_FALSE);CHKERRQ(ierr); 280f2cacb80SMatthew G. Knepley ierr = VecNorm(r, NORM_2, &res);CHKERRQ(ierr); 281f2cacb80SMatthew G. Knepley if (tol >= 0.0) { 282f2cacb80SMatthew G. Knepley if (res > tol) SETERRQ2(comm, PETSC_ERR_ARG_WRONG, "L_2 Residual %g exceeds tolerance %g", (double) res, (double) tol); 283f2cacb80SMatthew G. Knepley } else if (residual) { 284f2cacb80SMatthew G. Knepley *residual = res; 285f2cacb80SMatthew G. Knepley } else { 286f2cacb80SMatthew G. Knepley ierr = PetscPrintf(comm, "L_2 Residual: %g\n", (double)res);CHKERRQ(ierr); 287f2cacb80SMatthew G. Knepley ierr = VecChop(r, 1.0e-10);CHKERRQ(ierr); 2880c364540SMatthew G. Knepley ierr = PetscObjectCompose((PetscObject) r, "__Vec_bc_zero__", (PetscObject) dm);CHKERRQ(ierr); 289f2cacb80SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) r, "Initial Residual");CHKERRQ(ierr); 290f2cacb80SMatthew G. Knepley ierr = PetscObjectSetOptionsPrefix((PetscObject)r,"res_");CHKERRQ(ierr); 291f2cacb80SMatthew G. Knepley ierr = VecViewFromOptions(r, NULL, "-vec_view");CHKERRQ(ierr); 2920c364540SMatthew G. Knepley ierr = PetscObjectCompose((PetscObject) r, "__Vec_bc_zero__", NULL);CHKERRQ(ierr); 293f2cacb80SMatthew G. Knepley } 294f2cacb80SMatthew G. Knepley ierr = VecDestroy(&r);CHKERRQ(ierr); 295f2cacb80SMatthew G. Knepley PetscFunctionReturn(0); 296f2cacb80SMatthew G. Knepley } 297f2cacb80SMatthew G. Knepley 298f2cacb80SMatthew G. Knepley /*@C 299f2cacb80SMatthew G. Knepley DMTSCheckJacobian - Check the Jacobian of the exact solution against the residual using the Taylor Test 300f2cacb80SMatthew G. Knepley 301f2cacb80SMatthew G. Knepley Input Parameters: 302f2cacb80SMatthew G. Knepley + ts - the TS object 303f2cacb80SMatthew G. Knepley . dm - the DM 304f2cacb80SMatthew G. Knepley . t - the time 305f2cacb80SMatthew G. Knepley . u - a DM vector 306f2cacb80SMatthew G. Knepley . u_t - a DM vector 307f2cacb80SMatthew G. Knepley - tol - A tolerance for the check, or -1 to print the results instead 308f2cacb80SMatthew G. Knepley 309f2cacb80SMatthew G. Knepley Output Parameters: 310f2cacb80SMatthew G. Knepley + isLinear - Flag indicaing that the function looks linear, or NULL 311f2cacb80SMatthew G. Knepley - convRate - The rate of convergence of the linear model, or NULL 312f2cacb80SMatthew G. Knepley 313f2cacb80SMatthew G. Knepley Level: developer 314f2cacb80SMatthew G. Knepley 315f2cacb80SMatthew G. Knepley .seealso: DNTSCheckFromOptions(), DMTSCheckResidual(), DNSNESCheckFromOptions(), DMSNESCheckDiscretization(), DMSNESCheckResidual() 316f2cacb80SMatthew G. Knepley @*/ 317f2cacb80SMatthew G. Knepley PetscErrorCode DMTSCheckJacobian(TS ts, DM dm, PetscReal t, Vec u, Vec u_t, PetscReal tol, PetscBool *isLinear, PetscReal *convRate) 318f2cacb80SMatthew G. Knepley { 319f2cacb80SMatthew G. Knepley MPI_Comm comm; 320f2cacb80SMatthew G. Knepley PetscDS ds; 321f2cacb80SMatthew G. Knepley Mat J, M; 322f2cacb80SMatthew G. Knepley MatNullSpace nullspace; 323f2cacb80SMatthew G. Knepley PetscReal dt, shift, slope, intercept; 324f2cacb80SMatthew G. Knepley PetscBool hasJac, hasPrec, isLin = PETSC_FALSE; 325f2cacb80SMatthew G. Knepley PetscErrorCode ierr; 326f2cacb80SMatthew G. Knepley 327f2cacb80SMatthew G. Knepley PetscFunctionBegin; 328f2cacb80SMatthew G. Knepley PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 329f2cacb80SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 2); 330f2cacb80SMatthew G. Knepley PetscValidHeaderSpecific(u, VEC_CLASSID, 3); 331f2cacb80SMatthew G. Knepley if (isLinear) PetscValidBoolPointer(isLinear, 5); 332f2cacb80SMatthew G. Knepley if (convRate) PetscValidRealPointer(convRate, 5); 333f2cacb80SMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject) ts, &comm);CHKERRQ(ierr); 334f2cacb80SMatthew G. Knepley ierr = DMComputeExactSolution(dm, t, u, u_t);CHKERRQ(ierr); 335f2cacb80SMatthew G. Knepley /* Create and view matrices */ 336f2cacb80SMatthew G. Knepley ierr = TSGetTimeStep(ts, &dt);CHKERRQ(ierr); 337f2cacb80SMatthew G. Knepley shift = 1.0/dt; 338f2cacb80SMatthew G. Knepley ierr = DMCreateMatrix(dm, &J);CHKERRQ(ierr); 339f2cacb80SMatthew G. Knepley ierr = DMGetDS(dm, &ds);CHKERRQ(ierr); 340f2cacb80SMatthew G. Knepley ierr = PetscDSHasJacobian(ds, &hasJac);CHKERRQ(ierr); 341f2cacb80SMatthew G. Knepley ierr = PetscDSHasJacobianPreconditioner(ds, &hasPrec);CHKERRQ(ierr); 342f2cacb80SMatthew G. Knepley if (hasJac && hasPrec) { 343f2cacb80SMatthew G. Knepley ierr = DMCreateMatrix(dm, &M);CHKERRQ(ierr); 344f2cacb80SMatthew G. Knepley ierr = TSComputeIJacobian(ts, t, u, u_t, shift, J, M, PETSC_FALSE);CHKERRQ(ierr); 345f2cacb80SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) M, "Preconditioning Matrix");CHKERRQ(ierr); 346f2cacb80SMatthew G. Knepley ierr = PetscObjectSetOptionsPrefix((PetscObject) M, "jacpre_");CHKERRQ(ierr); 347f2cacb80SMatthew G. Knepley ierr = MatViewFromOptions(M, NULL, "-mat_view");CHKERRQ(ierr); 348f2cacb80SMatthew G. Knepley ierr = MatDestroy(&M);CHKERRQ(ierr); 349f2cacb80SMatthew G. Knepley } else { 350f2cacb80SMatthew G. Knepley ierr = TSComputeIJacobian(ts, t, u, u_t, shift, J, J, PETSC_FALSE);CHKERRQ(ierr); 351f2cacb80SMatthew G. Knepley } 352f2cacb80SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) J, "Jacobian");CHKERRQ(ierr); 353f2cacb80SMatthew G. Knepley ierr = PetscObjectSetOptionsPrefix((PetscObject) J, "jac_");CHKERRQ(ierr); 354f2cacb80SMatthew G. Knepley ierr = MatViewFromOptions(J, NULL, "-mat_view");CHKERRQ(ierr); 355f2cacb80SMatthew G. Knepley /* Check nullspace */ 356f2cacb80SMatthew G. Knepley ierr = MatGetNullSpace(J, &nullspace);CHKERRQ(ierr); 357f2cacb80SMatthew G. Knepley if (nullspace) { 358f2cacb80SMatthew G. Knepley PetscBool isNull; 359f2cacb80SMatthew G. Knepley ierr = MatNullSpaceTest(nullspace, J, &isNull);CHKERRQ(ierr); 360f2cacb80SMatthew G. Knepley if (!isNull) SETERRQ(comm, PETSC_ERR_PLIB, "The null space calculated for the system operator is invalid."); 361f2cacb80SMatthew G. Knepley } 362f2cacb80SMatthew G. Knepley /* Taylor test */ 363f2cacb80SMatthew G. Knepley { 364f2cacb80SMatthew G. Knepley PetscRandom rand; 365f2cacb80SMatthew G. Knepley Vec du, uhat, uhat_t, r, rhat, df; 366f2cacb80SMatthew G. Knepley PetscReal h; 367f2cacb80SMatthew G. Knepley PetscReal *es, *hs, *errors; 368f2cacb80SMatthew G. Knepley PetscReal hMax = 1.0, hMin = 1e-6, hMult = 0.1; 369f2cacb80SMatthew G. Knepley PetscInt Nv, v; 370f2cacb80SMatthew G. Knepley 371f2cacb80SMatthew G. Knepley /* Choose a perturbation direction */ 372f2cacb80SMatthew G. Knepley ierr = PetscRandomCreate(comm, &rand);CHKERRQ(ierr); 373f2cacb80SMatthew G. Knepley ierr = VecDuplicate(u, &du);CHKERRQ(ierr); 374f2cacb80SMatthew G. Knepley ierr = VecSetRandom(du, rand);CHKERRQ(ierr); 375f2cacb80SMatthew G. Knepley ierr = PetscRandomDestroy(&rand);CHKERRQ(ierr); 376f2cacb80SMatthew G. Knepley ierr = VecDuplicate(u, &df);CHKERRQ(ierr); 377f2cacb80SMatthew G. Knepley ierr = MatMult(J, du, df);CHKERRQ(ierr); 378f2cacb80SMatthew G. Knepley /* Evaluate residual at u, F(u), save in vector r */ 379f2cacb80SMatthew G. Knepley ierr = VecDuplicate(u, &r);CHKERRQ(ierr); 380f2cacb80SMatthew G. Knepley ierr = TSComputeIFunction(ts, t, u, u_t, r, PETSC_FALSE);CHKERRQ(ierr); 381f2cacb80SMatthew G. Knepley /* Look at the convergence of our Taylor approximation as we approach u */ 382f2cacb80SMatthew G. Knepley for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv); 383f2cacb80SMatthew G. Knepley ierr = PetscCalloc3(Nv, &es, Nv, &hs, Nv, &errors);CHKERRQ(ierr); 384f2cacb80SMatthew G. Knepley ierr = VecDuplicate(u, &uhat);CHKERRQ(ierr); 385f2cacb80SMatthew G. Knepley ierr = VecDuplicate(u, &uhat_t);CHKERRQ(ierr); 386f2cacb80SMatthew G. Knepley ierr = VecDuplicate(u, &rhat);CHKERRQ(ierr); 387f2cacb80SMatthew G. Knepley for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv) { 388f2cacb80SMatthew G. Knepley ierr = VecWAXPY(uhat, h, du, u);CHKERRQ(ierr); 389f2cacb80SMatthew G. Knepley ierr = VecWAXPY(uhat_t, h*shift, du, u_t);CHKERRQ(ierr); 390f2cacb80SMatthew 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 */ 391f2cacb80SMatthew G. Knepley ierr = TSComputeIFunction(ts, t, uhat, uhat_t, rhat, PETSC_FALSE);CHKERRQ(ierr); 392f2cacb80SMatthew G. Knepley ierr = VecAXPBYPCZ(rhat, -1.0, -h, 1.0, r, df);CHKERRQ(ierr); 393f2cacb80SMatthew G. Knepley ierr = VecNorm(rhat, NORM_2, &errors[Nv]);CHKERRQ(ierr); 394f2cacb80SMatthew G. Knepley 395f2cacb80SMatthew G. Knepley es[Nv] = PetscLog10Real(errors[Nv]); 396f2cacb80SMatthew G. Knepley hs[Nv] = PetscLog10Real(h); 397f2cacb80SMatthew G. Knepley } 398f2cacb80SMatthew G. Knepley ierr = VecDestroy(&uhat);CHKERRQ(ierr); 399f2cacb80SMatthew G. Knepley ierr = VecDestroy(&uhat_t);CHKERRQ(ierr); 400f2cacb80SMatthew G. Knepley ierr = VecDestroy(&rhat);CHKERRQ(ierr); 401f2cacb80SMatthew G. Knepley ierr = VecDestroy(&df);CHKERRQ(ierr); 402f2cacb80SMatthew G. Knepley ierr = VecDestroy(&r);CHKERRQ(ierr); 403f2cacb80SMatthew G. Knepley ierr = VecDestroy(&du);CHKERRQ(ierr); 404f2cacb80SMatthew G. Knepley for (v = 0; v < Nv; ++v) { 405f2cacb80SMatthew G. Knepley if ((tol >= 0) && (errors[v] > tol)) break; 406f2cacb80SMatthew G. Knepley else if (errors[v] > PETSC_SMALL) break; 407f2cacb80SMatthew G. Knepley } 408f2cacb80SMatthew G. Knepley if (v == Nv) isLin = PETSC_TRUE; 409f2cacb80SMatthew G. Knepley ierr = PetscLinearRegression(Nv, hs, es, &slope, &intercept);CHKERRQ(ierr); 410f2cacb80SMatthew G. Knepley ierr = PetscFree3(es, hs, errors);CHKERRQ(ierr); 411f2cacb80SMatthew G. Knepley /* Slope should be about 2 */ 412f2cacb80SMatthew G. Knepley if (tol >= 0) { 413f2cacb80SMatthew 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); 414f2cacb80SMatthew G. Knepley } else if (isLinear || convRate) { 415f2cacb80SMatthew G. Knepley if (isLinear) *isLinear = isLin; 416f2cacb80SMatthew G. Knepley if (convRate) *convRate = slope; 417f2cacb80SMatthew G. Knepley } else { 418f2cacb80SMatthew G. Knepley if (!isLin) {ierr = PetscPrintf(comm, "Taylor approximation converging at order %3.2f\n", (double) slope);CHKERRQ(ierr);} 419f2cacb80SMatthew G. Knepley else {ierr = PetscPrintf(comm, "Function appears to be linear\n");CHKERRQ(ierr);} 420f2cacb80SMatthew G. Knepley } 421f2cacb80SMatthew G. Knepley } 422f2cacb80SMatthew G. Knepley ierr = MatDestroy(&J);CHKERRQ(ierr); 423f2cacb80SMatthew G. Knepley PetscFunctionReturn(0); 424f2cacb80SMatthew G. Knepley } 425f2cacb80SMatthew G. Knepley 426f2cacb80SMatthew G. Knepley /*@C 427bee9a294SMatthew G. Knepley DMTSCheckFromOptions - Check the residual and Jacobian functions using the exact solution by outputting some diagnostic information 428bee9a294SMatthew G. Knepley 429bee9a294SMatthew G. Knepley Input Parameters: 430bee9a294SMatthew G. Knepley + ts - the TS object 4317f96f943SMatthew G. Knepley - u - representative TS vector 4327f96f943SMatthew G. Knepley 4337f96f943SMatthew G. Knepley Note: The user must call PetscDSSetExactSolution() beforehand 434bee9a294SMatthew G. Knepley 435bee9a294SMatthew G. Knepley Level: developer 436bee9a294SMatthew G. Knepley @*/ 4377f96f943SMatthew G. Knepley PetscErrorCode DMTSCheckFromOptions(TS ts, Vec u) 4387cdb2a12SMatthew G. Knepley { 4397cdb2a12SMatthew G. Knepley DM dm; 4407cdb2a12SMatthew G. Knepley SNES snes; 441f2cacb80SMatthew G. Knepley Vec sol, u_t; 442f2cacb80SMatthew G. Knepley PetscReal t; 4437cdb2a12SMatthew G. Knepley PetscBool check; 4447cdb2a12SMatthew G. Knepley PetscErrorCode ierr; 4457cdb2a12SMatthew G. Knepley 4467cdb2a12SMatthew G. Knepley PetscFunctionBegin; 447c5929fdfSBarry Smith ierr = PetscOptionsHasName(((PetscObject)ts)->options,((PetscObject)ts)->prefix, "-dmts_check", &check);CHKERRQ(ierr); 4487cdb2a12SMatthew G. Knepley if (!check) PetscFunctionReturn(0); 44972fd7fc8SMatthew G. Knepley ierr = VecDuplicate(u, &sol);CHKERRQ(ierr); 4506fb91e28SMatthew G. Knepley ierr = VecCopy(u, sol);CHKERRQ(ierr); 4516fb91e28SMatthew G. Knepley ierr = TSSetSolution(ts, u);CHKERRQ(ierr); 4527cdb2a12SMatthew G. Knepley ierr = TSGetDM(ts, &dm);CHKERRQ(ierr); 4537cdb2a12SMatthew G. Knepley ierr = TSSetUp(ts);CHKERRQ(ierr); 4547cdb2a12SMatthew G. Knepley ierr = TSGetSNES(ts, &snes);CHKERRQ(ierr); 4556fb91e28SMatthew G. Knepley ierr = SNESSetSolution(snes, u);CHKERRQ(ierr); 456f2cacb80SMatthew G. Knepley 457f2cacb80SMatthew G. Knepley ierr = TSGetTime(ts, &t);CHKERRQ(ierr); 458f2cacb80SMatthew G. Knepley ierr = DMSNESCheckDiscretization(snes, dm, t, sol, -1.0, NULL);CHKERRQ(ierr); 459f2cacb80SMatthew G. Knepley ierr = DMGetGlobalVector(dm, &u_t);CHKERRQ(ierr); 460f2cacb80SMatthew G. Knepley ierr = DMTSCheckResidual(ts, dm, t, sol, u_t, -1.0, NULL);CHKERRQ(ierr); 461f2cacb80SMatthew G. Knepley ierr = DMTSCheckJacobian(ts, dm, t, sol, u_t, -1.0, NULL, NULL);CHKERRQ(ierr); 462f2cacb80SMatthew G. Knepley ierr = DMRestoreGlobalVector(dm, &u_t);CHKERRQ(ierr); 463f2cacb80SMatthew G. Knepley 46472fd7fc8SMatthew G. Knepley ierr = VecDestroy(&sol);CHKERRQ(ierr); 4657cdb2a12SMatthew G. Knepley PetscFunctionReturn(0); 4667cdb2a12SMatthew G. Knepley } 467