108da532bSDmitry Karpeev #include <petscsnes.h> 2b45d2f2cSJed Brown #include <petsc-private/dmimpl.h> /*I "petscdm.h" I*/ 347c6ae99SBarry Smith 4732e2eb9SMatthew G Knepley PetscClassId DM_CLASSID; 567a56275SMatthew G Knepley PetscLogEvent DM_Convert, DM_GlobalToLocal, DM_LocalToGlobal; 667a56275SMatthew G Knepley 747c6ae99SBarry Smith #undef __FUNCT__ 8a4121054SBarry Smith #define __FUNCT__ "DMCreate" 9a4121054SBarry Smith /*@ 10de043629SMatthew G Knepley DMCreate - Creates an empty DM object. The type can then be set with DMSetType(). 11a4121054SBarry Smith 12a4121054SBarry Smith If you never call DMSetType() it will generate an 13a4121054SBarry Smith error when you try to use the vector. 14a4121054SBarry Smith 15a4121054SBarry Smith Collective on MPI_Comm 16a4121054SBarry Smith 17a4121054SBarry Smith Input Parameter: 18a4121054SBarry Smith . comm - The communicator for the DM object 19a4121054SBarry Smith 20a4121054SBarry Smith Output Parameter: 21a4121054SBarry Smith . dm - The DM object 22a4121054SBarry Smith 23a4121054SBarry Smith Level: beginner 24a4121054SBarry Smith 25a4121054SBarry Smith .seealso: DMSetType(), DMDA, DMSLICED, DMCOMPOSITE 26a4121054SBarry Smith @*/ 277087cfbeSBarry Smith PetscErrorCode DMCreate(MPI_Comm comm,DM *dm) 28a4121054SBarry Smith { 29a4121054SBarry Smith DM v; 30a4121054SBarry Smith PetscErrorCode ierr; 31a4121054SBarry Smith 32a4121054SBarry Smith PetscFunctionBegin; 331411c6eeSJed Brown PetscValidPointer(dm,2); 341411c6eeSJed Brown *dm = PETSC_NULL; 35a4121054SBarry Smith #ifndef PETSC_USE_DYNAMIC_LIBRARIES 36b84caa0eSBarry Smith ierr = VecInitializePackage(PETSC_NULL);CHKERRQ(ierr); 37b84caa0eSBarry Smith ierr = MatInitializePackage(PETSC_NULL);CHKERRQ(ierr); 38a4121054SBarry Smith ierr = DMInitializePackage(PETSC_NULL);CHKERRQ(ierr); 39a4121054SBarry Smith #endif 40a4121054SBarry Smith 413194b578SJed Brown ierr = PetscHeaderCreate(v, _p_DM, struct _DMOps, DM_CLASSID, -1, "DM", "Distribution Manager", "DM", comm, DMDestroy, DMView);CHKERRQ(ierr); 42a4121054SBarry Smith ierr = PetscMemzero(v->ops, sizeof(struct _DMOps));CHKERRQ(ierr); 431411c6eeSJed Brown 44e7c4fc90SDmitry Karpeev 45a89ea682SMatthew G Knepley v->workSize = 0; 46a89ea682SMatthew G Knepley v->workArray = PETSC_NULL; 471411c6eeSJed Brown v->ltogmap = PETSC_NULL; 481411c6eeSJed Brown v->ltogmapb = PETSC_NULL; 491411c6eeSJed Brown v->bs = 1; 50171400e9SBarry Smith v->coloringtype = IS_COLORING_GLOBAL; 511411c6eeSJed Brown 521411c6eeSJed Brown *dm = v; 53a4121054SBarry Smith PetscFunctionReturn(0); 54a4121054SBarry Smith } 55a4121054SBarry Smith 56a4121054SBarry Smith 57a4121054SBarry Smith #undef __FUNCT__ 589a42bb27SBarry Smith #define __FUNCT__ "DMSetVecType" 599a42bb27SBarry Smith /*@C 60564755cdSBarry Smith DMSetVecType - Sets the type of vector created with DMCreateLocalVector() and DMCreateGlobalVector() 619a42bb27SBarry Smith 62aa219208SBarry Smith Logically Collective on DMDA 639a42bb27SBarry Smith 649a42bb27SBarry Smith Input Parameter: 659a42bb27SBarry Smith + da - initial distributed array 668154be41SBarry Smith . ctype - the vector type, currently either VECSTANDARD or VECCUSP 679a42bb27SBarry Smith 689a42bb27SBarry Smith Options Database: 69dd85299cSBarry Smith . -dm_vec_type ctype 709a42bb27SBarry Smith 719a42bb27SBarry Smith Level: intermediate 729a42bb27SBarry Smith 73aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDA, DMDAInterpolationType, VecType 749a42bb27SBarry Smith @*/ 757087cfbeSBarry Smith PetscErrorCode DMSetVecType(DM da,const VecType ctype) 769a42bb27SBarry Smith { 779a42bb27SBarry Smith PetscErrorCode ierr; 789a42bb27SBarry Smith 799a42bb27SBarry Smith PetscFunctionBegin; 809a42bb27SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 819a42bb27SBarry Smith ierr = PetscFree(da->vectype);CHKERRQ(ierr); 829a42bb27SBarry Smith ierr = PetscStrallocpy(ctype,&da->vectype);CHKERRQ(ierr); 839a42bb27SBarry Smith PetscFunctionReturn(0); 849a42bb27SBarry Smith } 859a42bb27SBarry Smith 869a42bb27SBarry Smith #undef __FUNCT__ 87521d9a4cSLisandro Dalcin #define __FUNCT__ "DMSetMatType" 88521d9a4cSLisandro Dalcin /*@C 89521d9a4cSLisandro Dalcin DMSetMatType - Sets the type of matrix created with DMCreateMatrix() 90521d9a4cSLisandro Dalcin 91521d9a4cSLisandro Dalcin Logically Collective on DM 92521d9a4cSLisandro Dalcin 93521d9a4cSLisandro Dalcin Input Parameter: 94521d9a4cSLisandro Dalcin + dm - the DM context 95521d9a4cSLisandro Dalcin . ctype - the matrix type 96521d9a4cSLisandro Dalcin 97521d9a4cSLisandro Dalcin Options Database: 98521d9a4cSLisandro Dalcin . -dm_mat_type ctype 99521d9a4cSLisandro Dalcin 100521d9a4cSLisandro Dalcin Level: intermediate 101521d9a4cSLisandro Dalcin 102521d9a4cSLisandro Dalcin .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMCreateMatrix(), DMSetMatrixPreallocateOnly(), MatType 103521d9a4cSLisandro Dalcin @*/ 104521d9a4cSLisandro Dalcin PetscErrorCode DMSetMatType(DM dm,const MatType ctype) 105521d9a4cSLisandro Dalcin { 106521d9a4cSLisandro Dalcin PetscErrorCode ierr; 107521d9a4cSLisandro Dalcin PetscFunctionBegin; 108521d9a4cSLisandro Dalcin PetscValidHeaderSpecific(dm,DM_CLASSID,1); 109521d9a4cSLisandro Dalcin ierr = PetscFree(dm->mattype);CHKERRQ(ierr); 110521d9a4cSLisandro Dalcin ierr = PetscStrallocpy(ctype,&dm->mattype);CHKERRQ(ierr); 111521d9a4cSLisandro Dalcin PetscFunctionReturn(0); 112521d9a4cSLisandro Dalcin } 113521d9a4cSLisandro Dalcin 114521d9a4cSLisandro Dalcin #undef __FUNCT__ 1159a42bb27SBarry Smith #define __FUNCT__ "DMSetOptionsPrefix" 1169a42bb27SBarry Smith /*@C 1179a42bb27SBarry Smith DMSetOptionsPrefix - Sets the prefix used for searching for all 118aa219208SBarry Smith DMDA options in the database. 1199a42bb27SBarry Smith 120aa219208SBarry Smith Logically Collective on DMDA 1219a42bb27SBarry Smith 1229a42bb27SBarry Smith Input Parameter: 123aa219208SBarry Smith + da - the DMDA context 1249a42bb27SBarry Smith - prefix - the prefix to prepend to all option names 1259a42bb27SBarry Smith 1269a42bb27SBarry Smith Notes: 1279a42bb27SBarry Smith A hyphen (-) must NOT be given at the beginning of the prefix name. 1289a42bb27SBarry Smith The first character of all runtime options is AUTOMATICALLY the hyphen. 1299a42bb27SBarry Smith 1309a42bb27SBarry Smith Level: advanced 1319a42bb27SBarry Smith 132aa219208SBarry Smith .keywords: DMDA, set, options, prefix, database 1339a42bb27SBarry Smith 1349a42bb27SBarry Smith .seealso: DMSetFromOptions() 1359a42bb27SBarry Smith @*/ 1367087cfbeSBarry Smith PetscErrorCode DMSetOptionsPrefix(DM dm,const char prefix[]) 1379a42bb27SBarry Smith { 1389a42bb27SBarry Smith PetscErrorCode ierr; 1399a42bb27SBarry Smith 1409a42bb27SBarry Smith PetscFunctionBegin; 1419a42bb27SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1429a42bb27SBarry Smith ierr = PetscObjectSetOptionsPrefix((PetscObject)dm,prefix);CHKERRQ(ierr); 1439a42bb27SBarry Smith PetscFunctionReturn(0); 1449a42bb27SBarry Smith } 1459a42bb27SBarry Smith 1469a42bb27SBarry Smith #undef __FUNCT__ 14747c6ae99SBarry Smith #define __FUNCT__ "DMDestroy" 14847c6ae99SBarry Smith /*@ 149aa219208SBarry Smith DMDestroy - Destroys a vector packer or DMDA. 15047c6ae99SBarry Smith 15147c6ae99SBarry Smith Collective on DM 15247c6ae99SBarry Smith 15347c6ae99SBarry Smith Input Parameter: 15447c6ae99SBarry Smith . dm - the DM object to destroy 15547c6ae99SBarry Smith 15647c6ae99SBarry Smith Level: developer 15747c6ae99SBarry Smith 158e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix() 15947c6ae99SBarry Smith 16047c6ae99SBarry Smith @*/ 161fcfd50ebSBarry Smith PetscErrorCode DMDestroy(DM *dm) 16247c6ae99SBarry Smith { 163732e2eb9SMatthew G Knepley PetscInt i, cnt = 0; 164b17ce1afSJed Brown DMCoarsenHookLink link,next; 165*dfe15315SJed Brown DMNamedVecLink nlink,nnext; 16647c6ae99SBarry Smith PetscErrorCode ierr; 16747c6ae99SBarry Smith 16847c6ae99SBarry Smith PetscFunctionBegin; 1696bf464f9SBarry Smith if (!*dm) PetscFunctionReturn(0); 1706bf464f9SBarry Smith PetscValidHeaderSpecific((*dm),DM_CLASSID,1); 17187e657c6SBarry Smith 17287e657c6SBarry Smith /* count all the circular references of DM and its contained Vecs */ 173732e2eb9SMatthew G Knepley for (i=0; i<DM_MAX_WORK_VECTORS; i++) { 1746bf464f9SBarry Smith if ((*dm)->localin[i]) {cnt++;} 1756bf464f9SBarry Smith if ((*dm)->globalin[i]) {cnt++;} 176732e2eb9SMatthew G Knepley } 177*dfe15315SJed Brown for (nlink=(*dm)->namedglobal; nlink; nlink=nlink->next) cnt++; 17871cd77b2SBarry Smith if ((*dm)->x) { 17971cd77b2SBarry Smith PetscObject obj; 18071cd77b2SBarry Smith ierr = PetscObjectQuery((PetscObject)(*dm)->x,"DM",&obj);CHKERRQ(ierr); 18171cd77b2SBarry Smith if (obj == (PetscObject)*dm) cnt++; 18271cd77b2SBarry Smith } 183732e2eb9SMatthew G Knepley 1846bf464f9SBarry Smith if (--((PetscObject)(*dm))->refct - cnt > 0) {*dm = 0; PetscFunctionReturn(0);} 185732e2eb9SMatthew G Knepley /* 186732e2eb9SMatthew G Knepley Need this test because the dm references the vectors that 187732e2eb9SMatthew G Knepley reference the dm, so destroying the dm calls destroy on the 188732e2eb9SMatthew G Knepley vectors that cause another destroy on the dm 189732e2eb9SMatthew G Knepley */ 1906bf464f9SBarry Smith if (((PetscObject)(*dm))->refct < 0) PetscFunctionReturn(0); 1916bf464f9SBarry Smith ((PetscObject) (*dm))->refct = 0; 192732e2eb9SMatthew G Knepley for (i=0; i<DM_MAX_WORK_VECTORS; i++) { 1936bf464f9SBarry Smith if ((*dm)->localout[i]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Destroying a DM that has a local vector obtained with DMGetLocalVector()"); 1946bf464f9SBarry Smith ierr = VecDestroy(&(*dm)->localin[i]);CHKERRQ(ierr); 195732e2eb9SMatthew G Knepley } 196*dfe15315SJed Brown for (nlink=(*dm)->namedglobal; nlink; nlink=nnext) { /* Destroy the named vectors */ 197*dfe15315SJed Brown nnext = nlink->next; 198*dfe15315SJed Brown if (nlink->status != DMVEC_STATUS_IN) SETERRQ1(((PetscObject)*dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"DM still has Vec named '%s' checked out",nlink->name); 199*dfe15315SJed Brown ierr = PetscFree(nlink->name);CHKERRQ(ierr); 200*dfe15315SJed Brown ierr = VecDestroy(&nlink->X);CHKERRQ(ierr); 201*dfe15315SJed Brown ierr = PetscFree(nlink);CHKERRQ(ierr); 202*dfe15315SJed Brown } 203*dfe15315SJed Brown (*dm)->namedglobal = PETSC_NULL; 2041a266240SBarry Smith 205b17ce1afSJed Brown /* Destroy the list of hooks */ 206b17ce1afSJed Brown for (link=(*dm)->coarsenhook; link; link=next) { 207b17ce1afSJed Brown next = link->next; 208b17ce1afSJed Brown ierr = PetscFree(link);CHKERRQ(ierr); 209b17ce1afSJed Brown } 210b17ce1afSJed Brown (*dm)->coarsenhook = PETSC_NULL; 211b17ce1afSJed Brown 2121a266240SBarry Smith if ((*dm)->ctx && (*dm)->ctxdestroy) { 2131a266240SBarry Smith ierr = (*(*dm)->ctxdestroy)(&(*dm)->ctx);CHKERRQ(ierr); 2141a266240SBarry Smith } 21587e657c6SBarry Smith ierr = VecDestroy(&(*dm)->x);CHKERRQ(ierr); 21671cd77b2SBarry Smith ierr = MatFDColoringDestroy(&(*dm)->fd);CHKERRQ(ierr); 2174dcab191SBarry Smith ierr = DMClearGlobalVectors(*dm);CHKERRQ(ierr); 2186bf464f9SBarry Smith ierr = ISLocalToGlobalMappingDestroy(&(*dm)->ltogmap);CHKERRQ(ierr); 2196bf464f9SBarry Smith ierr = ISLocalToGlobalMappingDestroy(&(*dm)->ltogmapb);CHKERRQ(ierr); 2206bf464f9SBarry Smith ierr = PetscFree((*dm)->vectype);CHKERRQ(ierr); 221073dac72SJed Brown ierr = PetscFree((*dm)->mattype);CHKERRQ(ierr); 222a89ea682SMatthew G Knepley ierr = PetscFree((*dm)->workArray);CHKERRQ(ierr); 223732e2eb9SMatthew G Knepley /* if memory was published with AMS then destroy it */ 2246bf464f9SBarry Smith ierr = PetscObjectDepublish(*dm);CHKERRQ(ierr); 225732e2eb9SMatthew G Knepley 2266bf464f9SBarry Smith ierr = (*(*dm)->ops->destroy)(*dm);CHKERRQ(ierr); 2276bf464f9SBarry Smith ierr = PetscFree((*dm)->data);CHKERRQ(ierr); 228732e2eb9SMatthew G Knepley ierr = PetscHeaderDestroy(dm);CHKERRQ(ierr); 22947c6ae99SBarry Smith PetscFunctionReturn(0); 23047c6ae99SBarry Smith } 23147c6ae99SBarry Smith 23247c6ae99SBarry Smith #undef __FUNCT__ 233d7bf68aeSBarry Smith #define __FUNCT__ "DMSetUp" 234d7bf68aeSBarry Smith /*@ 235d7bf68aeSBarry Smith DMSetUp - sets up the data structures inside a DM object 236d7bf68aeSBarry Smith 237d7bf68aeSBarry Smith Collective on DM 238d7bf68aeSBarry Smith 239d7bf68aeSBarry Smith Input Parameter: 240d7bf68aeSBarry Smith . dm - the DM object to setup 241d7bf68aeSBarry Smith 242d7bf68aeSBarry Smith Level: developer 243d7bf68aeSBarry Smith 244e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix() 245d7bf68aeSBarry Smith 246d7bf68aeSBarry Smith @*/ 2477087cfbeSBarry Smith PetscErrorCode DMSetUp(DM dm) 248d7bf68aeSBarry Smith { 249d7bf68aeSBarry Smith PetscErrorCode ierr; 250d7bf68aeSBarry Smith 251d7bf68aeSBarry Smith PetscFunctionBegin; 252171400e9SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2538387afaaSJed Brown if (dm->setupcalled) PetscFunctionReturn(0); 254d7bf68aeSBarry Smith if (dm->ops->setup) { 255d7bf68aeSBarry Smith ierr = (*dm->ops->setup)(dm);CHKERRQ(ierr); 256d7bf68aeSBarry Smith } 2578387afaaSJed Brown dm->setupcalled = PETSC_TRUE; 258d7bf68aeSBarry Smith PetscFunctionReturn(0); 259d7bf68aeSBarry Smith } 260d7bf68aeSBarry Smith 261d7bf68aeSBarry Smith #undef __FUNCT__ 262d7bf68aeSBarry Smith #define __FUNCT__ "DMSetFromOptions" 263d7bf68aeSBarry Smith /*@ 264d7bf68aeSBarry Smith DMSetFromOptions - sets parameters in a DM from the options database 265d7bf68aeSBarry Smith 266d7bf68aeSBarry Smith Collective on DM 267d7bf68aeSBarry Smith 268d7bf68aeSBarry Smith Input Parameter: 269d7bf68aeSBarry Smith . dm - the DM object to set options for 270d7bf68aeSBarry Smith 271732e2eb9SMatthew G Knepley Options Database: 272dd85299cSBarry Smith + -dm_preallocate_only: Only preallocate the matrix for DMCreateMatrix(), but do not fill it with zeros 273dd85299cSBarry Smith . -dm_vec_type <type> type of vector to create inside DM 274171400e9SBarry Smith . -dm_mat_type <type> type of matrix to create inside DM 275171400e9SBarry Smith - -dm_coloring_type <global or ghosted> 276732e2eb9SMatthew G Knepley 277d7bf68aeSBarry Smith Level: developer 278d7bf68aeSBarry Smith 279e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix() 280d7bf68aeSBarry Smith 281d7bf68aeSBarry Smith @*/ 2827087cfbeSBarry Smith PetscErrorCode DMSetFromOptions(DM dm) 283d7bf68aeSBarry Smith { 28467ad5babSMatthew G Knepley PetscBool flg1 = PETSC_FALSE,flg2 = PETSC_FALSE,flg3 = PETSC_FALSE,flg4 = PETSC_FALSE,flg; 285d7bf68aeSBarry Smith PetscErrorCode ierr; 286f9ba7244SBarry Smith char typeName[256] = MATAIJ; 287d7bf68aeSBarry Smith 288d7bf68aeSBarry Smith PetscFunctionBegin; 289171400e9SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2903194b578SJed Brown ierr = PetscObjectOptionsBegin((PetscObject)dm);CHKERRQ(ierr); 29182fcb398SMatthew G Knepley ierr = PetscOptionsBool("-dm_view", "Information on DM", "DMView", flg1, &flg1, PETSC_NULL);CHKERRQ(ierr); 292c4721b0eSMatthew G Knepley ierr = PetscOptionsBool("-dm_view_detail", "Exhaustive mesh description", "DMView", flg2, &flg2, PETSC_NULL);CHKERRQ(ierr); 293c4721b0eSMatthew G Knepley ierr = PetscOptionsBool("-dm_view_vtk", "Output mesh in VTK format", "DMView", flg3, &flg3, PETSC_NULL);CHKERRQ(ierr); 29467ad5babSMatthew G Knepley ierr = PetscOptionsBool("-dm_view_latex", "Output mesh in LaTeX TikZ format", "DMView", flg4, &flg4, PETSC_NULL);CHKERRQ(ierr); 295073dac72SJed Brown ierr = PetscOptionsBool("-dm_preallocate_only","only preallocate matrix, but do not set column indices","DMSetMatrixPreallocateOnly",dm->prealloc_only,&dm->prealloc_only,PETSC_NULL);CHKERRQ(ierr); 296f9ba7244SBarry Smith ierr = PetscOptionsList("-dm_vec_type","Vector type used for created vectors","DMSetVecType",VecList,dm->vectype,typeName,256,&flg);CHKERRQ(ierr); 297f9ba7244SBarry Smith if (flg) { 298f9ba7244SBarry Smith ierr = DMSetVecType(dm,typeName);CHKERRQ(ierr); 299f9ba7244SBarry Smith } 300521d9a4cSLisandro Dalcin ierr = PetscOptionsList("-dm_mat_type","Matrix type used for created matrices","DMSetMatType",MatList,dm->mattype?dm->mattype:typeName,typeName,sizeof typeName,&flg);CHKERRQ(ierr); 301073dac72SJed Brown if (flg) { 302521d9a4cSLisandro Dalcin ierr = DMSetMatType(dm,typeName);CHKERRQ(ierr); 303073dac72SJed Brown } 3041b89239cSHong Zhang ierr = PetscOptionsEnum("-dm_is_coloring_type","Global or local coloring of Jacobian","ISColoringType",ISColoringTypes,(PetscEnum)dm->coloringtype,(PetscEnum*)&dm->coloringtype,PETSC_NULL);CHKERRQ(ierr); 305f9ba7244SBarry Smith if (dm->ops->setfromoptions) { 306f9ba7244SBarry Smith ierr = (*dm->ops->setfromoptions)(dm);CHKERRQ(ierr); 307f9ba7244SBarry Smith } 308f9ba7244SBarry Smith /* process any options handlers added with PetscObjectAddOptionsHandler() */ 309f9ba7244SBarry Smith ierr = PetscObjectProcessOptionsHandlers((PetscObject) dm);CHKERRQ(ierr); 31082fcb398SMatthew G Knepley ierr = PetscOptionsEnd();CHKERRQ(ierr); 31182fcb398SMatthew G Knepley if (flg1) { 31282fcb398SMatthew G Knepley ierr = DMView(dm, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 31382fcb398SMatthew G Knepley } 314c4721b0eSMatthew G Knepley if (flg2) { 315c4721b0eSMatthew G Knepley PetscViewer viewer; 316c4721b0eSMatthew G Knepley 317c4721b0eSMatthew G Knepley ierr = PetscViewerCreate(((PetscObject) dm)->comm, &viewer);CHKERRQ(ierr); 318c4721b0eSMatthew G Knepley ierr = PetscViewerSetType(viewer, PETSCVIEWERASCII);CHKERRQ(ierr); 319c4721b0eSMatthew G Knepley ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_INFO_DETAIL);CHKERRQ(ierr); 320c4721b0eSMatthew G Knepley ierr = DMView(dm, viewer);CHKERRQ(ierr); 321c4721b0eSMatthew G Knepley ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 322c4721b0eSMatthew G Knepley } 323c4721b0eSMatthew G Knepley if (flg3) { 324c4721b0eSMatthew G Knepley PetscViewer viewer; 325c4721b0eSMatthew G Knepley 326c4721b0eSMatthew G Knepley ierr = PetscViewerCreate(((PetscObject) dm)->comm, &viewer);CHKERRQ(ierr); 327c4721b0eSMatthew G Knepley ierr = PetscViewerSetType(viewer, PETSCVIEWERASCII);CHKERRQ(ierr); 328c4721b0eSMatthew G Knepley ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_VTK);CHKERRQ(ierr); 329c4721b0eSMatthew G Knepley ierr = PetscViewerFileSetName(viewer, "mesh.vtk");CHKERRQ(ierr); 330c4721b0eSMatthew G Knepley ierr = DMView(dm, viewer);CHKERRQ(ierr); 331c4721b0eSMatthew G Knepley ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 332c4721b0eSMatthew G Knepley } 33367ad5babSMatthew G Knepley if (flg4) { 33467ad5babSMatthew G Knepley PetscViewer viewer; 33567ad5babSMatthew G Knepley 33667ad5babSMatthew G Knepley ierr = PetscViewerCreate(((PetscObject) dm)->comm, &viewer);CHKERRQ(ierr); 33767ad5babSMatthew G Knepley ierr = PetscViewerSetType(viewer, PETSCVIEWERASCII);CHKERRQ(ierr); 33867ad5babSMatthew G Knepley ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_LATEX);CHKERRQ(ierr); 33967ad5babSMatthew G Knepley ierr = PetscViewerFileSetName(viewer, "mesh.tex");CHKERRQ(ierr); 34067ad5babSMatthew G Knepley ierr = DMView(dm, viewer);CHKERRQ(ierr); 34167ad5babSMatthew G Knepley ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 34267ad5babSMatthew G Knepley } 343d7bf68aeSBarry Smith PetscFunctionReturn(0); 344d7bf68aeSBarry Smith } 345d7bf68aeSBarry Smith 346d7bf68aeSBarry Smith #undef __FUNCT__ 34747c6ae99SBarry Smith #define __FUNCT__ "DMView" 348fc9bc008SSatish Balay /*@C 349aa219208SBarry Smith DMView - Views a vector packer or DMDA. 35047c6ae99SBarry Smith 35147c6ae99SBarry Smith Collective on DM 35247c6ae99SBarry Smith 35347c6ae99SBarry Smith Input Parameter: 35447c6ae99SBarry Smith + dm - the DM object to view 35547c6ae99SBarry Smith - v - the viewer 35647c6ae99SBarry Smith 35747c6ae99SBarry Smith Level: developer 35847c6ae99SBarry Smith 359e727c939SJed Brown .seealso DMDestroy(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix() 36047c6ae99SBarry Smith 36147c6ae99SBarry Smith @*/ 3627087cfbeSBarry Smith PetscErrorCode DMView(DM dm,PetscViewer v) 36347c6ae99SBarry Smith { 36447c6ae99SBarry Smith PetscErrorCode ierr; 36547c6ae99SBarry Smith 36647c6ae99SBarry Smith PetscFunctionBegin; 367171400e9SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 3683014e516SBarry Smith if (!v) { 3693014e516SBarry Smith ierr = PetscViewerASCIIGetStdout(((PetscObject)dm)->comm,&v);CHKERRQ(ierr); 3703014e516SBarry Smith } 3710c010503SBarry Smith if (dm->ops->view) { 3720c010503SBarry Smith ierr = (*dm->ops->view)(dm,v);CHKERRQ(ierr); 37347c6ae99SBarry Smith } 37447c6ae99SBarry Smith PetscFunctionReturn(0); 37547c6ae99SBarry Smith } 37647c6ae99SBarry Smith 37747c6ae99SBarry Smith #undef __FUNCT__ 37847c6ae99SBarry Smith #define __FUNCT__ "DMCreateGlobalVector" 37947c6ae99SBarry Smith /*@ 380aa219208SBarry Smith DMCreateGlobalVector - Creates a global vector from a DMDA or DMComposite object 38147c6ae99SBarry Smith 38247c6ae99SBarry Smith Collective on DM 38347c6ae99SBarry Smith 38447c6ae99SBarry Smith Input Parameter: 38547c6ae99SBarry Smith . dm - the DM object 38647c6ae99SBarry Smith 38747c6ae99SBarry Smith Output Parameter: 38847c6ae99SBarry Smith . vec - the global vector 38947c6ae99SBarry Smith 390073dac72SJed Brown Level: beginner 39147c6ae99SBarry Smith 392e727c939SJed Brown .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix() 39347c6ae99SBarry Smith 39447c6ae99SBarry Smith @*/ 3957087cfbeSBarry Smith PetscErrorCode DMCreateGlobalVector(DM dm,Vec *vec) 39647c6ae99SBarry Smith { 39747c6ae99SBarry Smith PetscErrorCode ierr; 39847c6ae99SBarry Smith 39947c6ae99SBarry Smith PetscFunctionBegin; 400171400e9SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 40147c6ae99SBarry Smith ierr = (*dm->ops->createglobalvector)(dm,vec);CHKERRQ(ierr); 40247c6ae99SBarry Smith PetscFunctionReturn(0); 40347c6ae99SBarry Smith } 40447c6ae99SBarry Smith 40547c6ae99SBarry Smith #undef __FUNCT__ 40647c6ae99SBarry Smith #define __FUNCT__ "DMCreateLocalVector" 40747c6ae99SBarry Smith /*@ 408aa219208SBarry Smith DMCreateLocalVector - Creates a local vector from a DMDA or DMComposite object 40947c6ae99SBarry Smith 41047c6ae99SBarry Smith Not Collective 41147c6ae99SBarry Smith 41247c6ae99SBarry Smith Input Parameter: 41347c6ae99SBarry Smith . dm - the DM object 41447c6ae99SBarry Smith 41547c6ae99SBarry Smith Output Parameter: 41647c6ae99SBarry Smith . vec - the local vector 41747c6ae99SBarry Smith 418073dac72SJed Brown Level: beginner 41947c6ae99SBarry Smith 420e727c939SJed Brown .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix() 42147c6ae99SBarry Smith 42247c6ae99SBarry Smith @*/ 4237087cfbeSBarry Smith PetscErrorCode DMCreateLocalVector(DM dm,Vec *vec) 42447c6ae99SBarry Smith { 42547c6ae99SBarry Smith PetscErrorCode ierr; 42647c6ae99SBarry Smith 42747c6ae99SBarry Smith PetscFunctionBegin; 428171400e9SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 42947c6ae99SBarry Smith ierr = (*dm->ops->createlocalvector)(dm,vec);CHKERRQ(ierr); 43047c6ae99SBarry Smith PetscFunctionReturn(0); 43147c6ae99SBarry Smith } 43247c6ae99SBarry Smith 43347c6ae99SBarry Smith #undef __FUNCT__ 4341411c6eeSJed Brown #define __FUNCT__ "DMGetLocalToGlobalMapping" 4351411c6eeSJed Brown /*@ 4361411c6eeSJed Brown DMGetLocalToGlobalMapping - Accesses the local-to-global mapping in a DM. 4371411c6eeSJed Brown 4381411c6eeSJed Brown Collective on DM 4391411c6eeSJed Brown 4401411c6eeSJed Brown Input Parameter: 4411411c6eeSJed Brown . dm - the DM that provides the mapping 4421411c6eeSJed Brown 4431411c6eeSJed Brown Output Parameter: 4441411c6eeSJed Brown . ltog - the mapping 4451411c6eeSJed Brown 4461411c6eeSJed Brown Level: intermediate 4471411c6eeSJed Brown 4481411c6eeSJed Brown Notes: 4491411c6eeSJed Brown This mapping can then be used by VecSetLocalToGlobalMapping() or 4501411c6eeSJed Brown MatSetLocalToGlobalMapping(). 4511411c6eeSJed Brown 4521411c6eeSJed Brown .seealso: DMCreateLocalVector(), DMGetLocalToGlobalMappingBlock() 4531411c6eeSJed Brown @*/ 4547087cfbeSBarry Smith PetscErrorCode DMGetLocalToGlobalMapping(DM dm,ISLocalToGlobalMapping *ltog) 4551411c6eeSJed Brown { 4561411c6eeSJed Brown PetscErrorCode ierr; 4571411c6eeSJed Brown 4581411c6eeSJed Brown PetscFunctionBegin; 4591411c6eeSJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 4601411c6eeSJed Brown PetscValidPointer(ltog,2); 4611411c6eeSJed Brown if (!dm->ltogmap) { 4621411c6eeSJed Brown if (!dm->ops->createlocaltoglobalmapping) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"DM can not create LocalToGlobalMapping"); 4631411c6eeSJed Brown ierr = (*dm->ops->createlocaltoglobalmapping)(dm);CHKERRQ(ierr); 4641411c6eeSJed Brown } 4651411c6eeSJed Brown *ltog = dm->ltogmap; 4661411c6eeSJed Brown PetscFunctionReturn(0); 4671411c6eeSJed Brown } 4681411c6eeSJed Brown 4691411c6eeSJed Brown #undef __FUNCT__ 4701411c6eeSJed Brown #define __FUNCT__ "DMGetLocalToGlobalMappingBlock" 4711411c6eeSJed Brown /*@ 4721411c6eeSJed Brown DMGetLocalToGlobalMappingBlock - Accesses the blocked local-to-global mapping in a DM. 4731411c6eeSJed Brown 4741411c6eeSJed Brown Collective on DM 4751411c6eeSJed Brown 4761411c6eeSJed Brown Input Parameter: 4771411c6eeSJed Brown . da - the distributed array that provides the mapping 4781411c6eeSJed Brown 4791411c6eeSJed Brown Output Parameter: 4801411c6eeSJed Brown . ltog - the block mapping 4811411c6eeSJed Brown 4821411c6eeSJed Brown Level: intermediate 4831411c6eeSJed Brown 4841411c6eeSJed Brown Notes: 4851411c6eeSJed Brown This mapping can then be used by VecSetLocalToGlobalMappingBlock() or 4861411c6eeSJed Brown MatSetLocalToGlobalMappingBlock(). 4871411c6eeSJed Brown 4881411c6eeSJed Brown .seealso: DMCreateLocalVector(), DMGetLocalToGlobalMapping(), DMGetBlockSize(), VecSetBlockSize(), MatSetBlockSize() 4891411c6eeSJed Brown @*/ 4907087cfbeSBarry Smith PetscErrorCode DMGetLocalToGlobalMappingBlock(DM dm,ISLocalToGlobalMapping *ltog) 4911411c6eeSJed Brown { 4921411c6eeSJed Brown PetscErrorCode ierr; 4931411c6eeSJed Brown 4941411c6eeSJed Brown PetscFunctionBegin; 4951411c6eeSJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 4961411c6eeSJed Brown PetscValidPointer(ltog,2); 4971411c6eeSJed Brown if (!dm->ltogmapb) { 4981411c6eeSJed Brown PetscInt bs; 4991411c6eeSJed Brown ierr = DMGetBlockSize(dm,&bs);CHKERRQ(ierr); 5001411c6eeSJed Brown if (bs > 1) { 5011411c6eeSJed Brown if (!dm->ops->createlocaltoglobalmappingblock) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"DM can not create LocalToGlobalMappingBlock"); 5021411c6eeSJed Brown ierr = (*dm->ops->createlocaltoglobalmappingblock)(dm);CHKERRQ(ierr); 5031411c6eeSJed Brown } else { 5041411c6eeSJed Brown ierr = DMGetLocalToGlobalMapping(dm,&dm->ltogmapb);CHKERRQ(ierr); 5051411c6eeSJed Brown ierr = PetscObjectReference((PetscObject)dm->ltogmapb);CHKERRQ(ierr); 5061411c6eeSJed Brown } 5071411c6eeSJed Brown } 5081411c6eeSJed Brown *ltog = dm->ltogmapb; 5091411c6eeSJed Brown PetscFunctionReturn(0); 5101411c6eeSJed Brown } 5111411c6eeSJed Brown 5121411c6eeSJed Brown #undef __FUNCT__ 5131411c6eeSJed Brown #define __FUNCT__ "DMGetBlockSize" 5141411c6eeSJed Brown /*@ 5151411c6eeSJed Brown DMGetBlockSize - Gets the inherent block size associated with a DM 5161411c6eeSJed Brown 5171411c6eeSJed Brown Not Collective 5181411c6eeSJed Brown 5191411c6eeSJed Brown Input Parameter: 5201411c6eeSJed Brown . dm - the DM with block structure 5211411c6eeSJed Brown 5221411c6eeSJed Brown Output Parameter: 5231411c6eeSJed Brown . bs - the block size, 1 implies no exploitable block structure 5241411c6eeSJed Brown 5251411c6eeSJed Brown Level: intermediate 5261411c6eeSJed Brown 5271411c6eeSJed Brown .seealso: ISCreateBlock(), VecSetBlockSize(), MatSetBlockSize(), DMGetLocalToGlobalMappingBlock() 5281411c6eeSJed Brown @*/ 5297087cfbeSBarry Smith PetscErrorCode DMGetBlockSize(DM dm,PetscInt *bs) 5301411c6eeSJed Brown { 5311411c6eeSJed Brown PetscFunctionBegin; 5321411c6eeSJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 5331411c6eeSJed Brown PetscValidPointer(bs,2); 5341411c6eeSJed Brown if (dm->bs < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"DM does not have enough information to provide a block size yet"); 5351411c6eeSJed Brown *bs = dm->bs; 5361411c6eeSJed Brown PetscFunctionReturn(0); 5371411c6eeSJed Brown } 5381411c6eeSJed Brown 5391411c6eeSJed Brown #undef __FUNCT__ 540e727c939SJed Brown #define __FUNCT__ "DMCreateInterpolation" 54147c6ae99SBarry Smith /*@ 542e727c939SJed Brown DMCreateInterpolation - Gets interpolation matrix between two DMDA or DMComposite objects 54347c6ae99SBarry Smith 54447c6ae99SBarry Smith Collective on DM 54547c6ae99SBarry Smith 54647c6ae99SBarry Smith Input Parameter: 54747c6ae99SBarry Smith + dm1 - the DM object 54847c6ae99SBarry Smith - dm2 - the second, finer DM object 54947c6ae99SBarry Smith 55047c6ae99SBarry Smith Output Parameter: 55147c6ae99SBarry Smith + mat - the interpolation 55247c6ae99SBarry Smith - vec - the scaling (optional) 55347c6ae99SBarry Smith 55447c6ae99SBarry Smith Level: developer 55547c6ae99SBarry Smith 55685afcc9aSBarry Smith Notes: For DMDA objects this only works for "uniform refinement", that is the refined mesh was obtained DMRefine() or the coarse mesh was obtained by 55785afcc9aSBarry Smith DMCoarsen(). The coordinates set into the DMDA are completely ignored in computing the interpolation. 558d52bd9f3SBarry Smith 559d52bd9f3SBarry Smith For DMDA objects you can use this interpolation (more precisely the interpolation from the DMDAGetCoordinateDA()) to interpolate the mesh coordinate vectors 560d52bd9f3SBarry Smith EXCEPT in the periodic case where it does not make sense since the coordinate vectors are not periodic. 56185afcc9aSBarry Smith 56285afcc9aSBarry Smith 563e727c939SJed Brown .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateColoring(), DMCreateMatrix(), DMRefine(), DMCoarsen() 56447c6ae99SBarry Smith 56547c6ae99SBarry Smith @*/ 566e727c939SJed Brown PetscErrorCode DMCreateInterpolation(DM dm1,DM dm2,Mat *mat,Vec *vec) 56747c6ae99SBarry Smith { 56847c6ae99SBarry Smith PetscErrorCode ierr; 56947c6ae99SBarry Smith 57047c6ae99SBarry Smith PetscFunctionBegin; 571171400e9SBarry Smith PetscValidHeaderSpecific(dm1,DM_CLASSID,1); 572171400e9SBarry Smith PetscValidHeaderSpecific(dm2,DM_CLASSID,2); 57325296bd5SBarry Smith ierr = (*dm1->ops->createinterpolation)(dm1,dm2,mat,vec);CHKERRQ(ierr); 57447c6ae99SBarry Smith PetscFunctionReturn(0); 57547c6ae99SBarry Smith } 57647c6ae99SBarry Smith 57747c6ae99SBarry Smith #undef __FUNCT__ 578e727c939SJed Brown #define __FUNCT__ "DMCreateInjection" 57947c6ae99SBarry Smith /*@ 580e727c939SJed Brown DMCreateInjection - Gets injection matrix between two DMDA or DMComposite objects 58147c6ae99SBarry Smith 58247c6ae99SBarry Smith Collective on DM 58347c6ae99SBarry Smith 58447c6ae99SBarry Smith Input Parameter: 58547c6ae99SBarry Smith + dm1 - the DM object 58647c6ae99SBarry Smith - dm2 - the second, finer DM object 58747c6ae99SBarry Smith 58847c6ae99SBarry Smith Output Parameter: 58947c6ae99SBarry Smith . ctx - the injection 59047c6ae99SBarry Smith 59147c6ae99SBarry Smith Level: developer 59247c6ae99SBarry Smith 59385afcc9aSBarry Smith Notes: For DMDA objects this only works for "uniform refinement", that is the refined mesh was obtained DMRefine() or the coarse mesh was obtained by 59485afcc9aSBarry Smith DMCoarsen(). The coordinates set into the DMDA are completely ignored in computing the injection. 59585afcc9aSBarry Smith 596e727c939SJed Brown .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateColoring(), DMCreateMatrix(), DMCreateInterpolation() 59747c6ae99SBarry Smith 59847c6ae99SBarry Smith @*/ 599e727c939SJed Brown PetscErrorCode DMCreateInjection(DM dm1,DM dm2,VecScatter *ctx) 60047c6ae99SBarry Smith { 60147c6ae99SBarry Smith PetscErrorCode ierr; 60247c6ae99SBarry Smith 60347c6ae99SBarry Smith PetscFunctionBegin; 604171400e9SBarry Smith PetscValidHeaderSpecific(dm1,DM_CLASSID,1); 605171400e9SBarry Smith PetscValidHeaderSpecific(dm2,DM_CLASSID,2); 60647c6ae99SBarry Smith ierr = (*dm1->ops->getinjection)(dm1,dm2,ctx);CHKERRQ(ierr); 60747c6ae99SBarry Smith PetscFunctionReturn(0); 60847c6ae99SBarry Smith } 60947c6ae99SBarry Smith 61047c6ae99SBarry Smith #undef __FUNCT__ 611e727c939SJed Brown #define __FUNCT__ "DMCreateColoring" 612d1e2c406SBarry Smith /*@C 613e727c939SJed Brown DMCreateColoring - Gets coloring for a DMDA or DMComposite 61447c6ae99SBarry Smith 61547c6ae99SBarry Smith Collective on DM 61647c6ae99SBarry Smith 61747c6ae99SBarry Smith Input Parameter: 61847c6ae99SBarry Smith + dm - the DM object 61947c6ae99SBarry Smith . ctype - IS_COLORING_GHOSTED or IS_COLORING_GLOBAL 62047c6ae99SBarry Smith - matype - either MATAIJ or MATBAIJ 62147c6ae99SBarry Smith 62247c6ae99SBarry Smith Output Parameter: 62347c6ae99SBarry Smith . coloring - the coloring 62447c6ae99SBarry Smith 62547c6ae99SBarry Smith Level: developer 62647c6ae99SBarry Smith 627e727c939SJed Brown .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateMatrix() 62847c6ae99SBarry Smith 629aab9d709SJed Brown @*/ 630e727c939SJed Brown PetscErrorCode DMCreateColoring(DM dm,ISColoringType ctype,const MatType mtype,ISColoring *coloring) 63147c6ae99SBarry Smith { 63247c6ae99SBarry Smith PetscErrorCode ierr; 63347c6ae99SBarry Smith 63447c6ae99SBarry Smith PetscFunctionBegin; 635171400e9SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 63647c6ae99SBarry Smith if (!dm->ops->getcoloring) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No coloring for this type of DM yet"); 63747c6ae99SBarry Smith ierr = (*dm->ops->getcoloring)(dm,ctype,mtype,coloring);CHKERRQ(ierr); 63847c6ae99SBarry Smith PetscFunctionReturn(0); 63947c6ae99SBarry Smith } 64047c6ae99SBarry Smith 64147c6ae99SBarry Smith #undef __FUNCT__ 642950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix" 64347c6ae99SBarry Smith /*@C 644950540a4SJed Brown DMCreateMatrix - Gets empty Jacobian for a DMDA or DMComposite 64547c6ae99SBarry Smith 64647c6ae99SBarry Smith Collective on DM 64747c6ae99SBarry Smith 64847c6ae99SBarry Smith Input Parameter: 64947c6ae99SBarry Smith + dm - the DM object 65047c6ae99SBarry Smith - mtype - Supported types are MATSEQAIJ, MATMPIAIJ, MATSEQBAIJ, MATMPIBAIJ, or 65194013140SBarry Smith any type which inherits from one of these (such as MATAIJ) 65247c6ae99SBarry Smith 65347c6ae99SBarry Smith Output Parameter: 65447c6ae99SBarry Smith . mat - the empty Jacobian 65547c6ae99SBarry Smith 656073dac72SJed Brown Level: beginner 65747c6ae99SBarry Smith 65894013140SBarry Smith Notes: This properly preallocates the number of nonzeros in the sparse matrix so you 65994013140SBarry Smith do not need to do it yourself. 66094013140SBarry Smith 66194013140SBarry Smith By default it also sets the nonzero structure and puts in the zero entries. To prevent setting 662aa219208SBarry Smith the nonzero pattern call DMDASetMatPreallocateOnly() 66394013140SBarry Smith 66494013140SBarry Smith For structured grid problems, when you call MatView() on this matrix it is displayed using the global natural ordering, NOT in the ordering used 66594013140SBarry Smith internally by PETSc. 66694013140SBarry Smith 66794013140SBarry Smith For structured grid problems, in general it is easiest to use MatSetValuesStencil() or MatSetValuesLocal() to put values into the matrix because MatSetValues() requires 668aa219208SBarry Smith the indices for the global numbering for DMDAs which is complicated. 66994013140SBarry Smith 670e727c939SJed Brown .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation() 67147c6ae99SBarry Smith 672aab9d709SJed Brown @*/ 673950540a4SJed Brown PetscErrorCode DMCreateMatrix(DM dm,const MatType mtype,Mat *mat) 67447c6ae99SBarry Smith { 67547c6ae99SBarry Smith PetscErrorCode ierr; 67647c6ae99SBarry Smith 67747c6ae99SBarry Smith PetscFunctionBegin; 678171400e9SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 679235683edSBarry Smith #ifndef PETSC_USE_DYNAMIC_LIBRARIES 680235683edSBarry Smith ierr = MatInitializePackage(PETSC_NULL);CHKERRQ(ierr); 681235683edSBarry Smith #endif 682c7b7c8a4SJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 683c7b7c8a4SJed Brown PetscValidPointer(mat,3); 684073dac72SJed Brown if (dm->mattype) { 68525296bd5SBarry Smith ierr = (*dm->ops->creatematrix)(dm,dm->mattype,mat);CHKERRQ(ierr); 686073dac72SJed Brown } else { 68725296bd5SBarry Smith ierr = (*dm->ops->creatematrix)(dm,mtype,mat);CHKERRQ(ierr); 688c7b7c8a4SJed Brown } 68947c6ae99SBarry Smith PetscFunctionReturn(0); 69047c6ae99SBarry Smith } 69147c6ae99SBarry Smith 69247c6ae99SBarry Smith #undef __FUNCT__ 693732e2eb9SMatthew G Knepley #define __FUNCT__ "DMSetMatrixPreallocateOnly" 694732e2eb9SMatthew G Knepley /*@ 695950540a4SJed Brown DMSetMatrixPreallocateOnly - When DMCreateMatrix() is called the matrix will be properly 696732e2eb9SMatthew G Knepley preallocated but the nonzero structure and zero values will not be set. 697732e2eb9SMatthew G Knepley 698732e2eb9SMatthew G Knepley Logically Collective on DMDA 699732e2eb9SMatthew G Knepley 700732e2eb9SMatthew G Knepley Input Parameter: 701732e2eb9SMatthew G Knepley + dm - the DM 702732e2eb9SMatthew G Knepley - only - PETSC_TRUE if only want preallocation 703732e2eb9SMatthew G Knepley 704732e2eb9SMatthew G Knepley Level: developer 705950540a4SJed Brown .seealso DMCreateMatrix() 706732e2eb9SMatthew G Knepley @*/ 707732e2eb9SMatthew G Knepley PetscErrorCode DMSetMatrixPreallocateOnly(DM dm, PetscBool only) 708732e2eb9SMatthew G Knepley { 709732e2eb9SMatthew G Knepley PetscFunctionBegin; 710732e2eb9SMatthew G Knepley PetscValidHeaderSpecific(dm,DM_CLASSID,1); 711732e2eb9SMatthew G Knepley dm->prealloc_only = only; 712732e2eb9SMatthew G Knepley PetscFunctionReturn(0); 713732e2eb9SMatthew G Knepley } 714732e2eb9SMatthew G Knepley 715732e2eb9SMatthew G Knepley #undef __FUNCT__ 716a89ea682SMatthew G Knepley #define __FUNCT__ "DMGetWorkArray" 717a89ea682SMatthew G Knepley /*@C 718a89ea682SMatthew G Knepley DMGetWorkArray - Gets a work array guaranteed to be at least the input size 719a89ea682SMatthew G Knepley 720a89ea682SMatthew G Knepley Not Collective 721a89ea682SMatthew G Knepley 722a89ea682SMatthew G Knepley Input Parameters: 723a89ea682SMatthew G Knepley + dm - the DM object 724a89ea682SMatthew G Knepley - size - The minium size 725a89ea682SMatthew G Knepley 726a89ea682SMatthew G Knepley Output Parameter: 727a89ea682SMatthew G Knepley . array - the work array 728a89ea682SMatthew G Knepley 729a89ea682SMatthew G Knepley Level: developer 730a89ea682SMatthew G Knepley 731a89ea682SMatthew G Knepley .seealso DMDestroy(), DMCreate() 732a89ea682SMatthew G Knepley @*/ 733a89ea682SMatthew G Knepley PetscErrorCode DMGetWorkArray(DM dm,PetscInt size,PetscScalar **array) 734a89ea682SMatthew G Knepley { 735a89ea682SMatthew G Knepley PetscErrorCode ierr; 736a89ea682SMatthew G Knepley 737a89ea682SMatthew G Knepley PetscFunctionBegin; 738a89ea682SMatthew G Knepley PetscValidHeaderSpecific(dm,DM_CLASSID,1); 739a89ea682SMatthew G Knepley PetscValidPointer(array,3); 740a89ea682SMatthew G Knepley if (size > dm->workSize) { 741a89ea682SMatthew G Knepley dm->workSize = size; 742a89ea682SMatthew G Knepley ierr = PetscFree(dm->workArray);CHKERRQ(ierr); 743a89ea682SMatthew G Knepley ierr = PetscMalloc(dm->workSize * sizeof(PetscScalar), &dm->workArray);CHKERRQ(ierr); 744a89ea682SMatthew G Knepley } 745a89ea682SMatthew G Knepley *array = dm->workArray; 746a89ea682SMatthew G Knepley PetscFunctionReturn(0); 747a89ea682SMatthew G Knepley } 748a89ea682SMatthew G Knepley 7494d343eeaSMatthew G Knepley #undef __FUNCT__ 750e7c4fc90SDmitry Karpeev #define __FUNCT__ "DMCreateDecompositionDM" 751e7c4fc90SDmitry Karpeev /*@C 752e7c4fc90SDmitry Karpeev DMCreateDecompositionDM - creates a DM that defines a decomposition of the original DM. 753e7c4fc90SDmitry Karpeev 754e7c4fc90SDmitry Karpeev Not Collective 755e7c4fc90SDmitry Karpeev 756e7c4fc90SDmitry Karpeev Input Parameters: 757e7c4fc90SDmitry Karpeev + dm - the DM object 758e7c4fc90SDmitry Karpeev - name - the name of the decomposition 759e7c4fc90SDmitry Karpeev 760e7c4fc90SDmitry Karpeev Output Parameter: 761e7c4fc90SDmitry Karpeev . ddm - the decomposition DM (PETSC_NULL, if no such decomposition is known) 762e7c4fc90SDmitry Karpeev 763e7c4fc90SDmitry Karpeev Level: advanced 764e7c4fc90SDmitry Karpeev 765e7c4fc90SDmitry Karpeev .seealso DMDestroy(), DMCreate(), DMCreateDecomposition() 766e7c4fc90SDmitry Karpeev @*/ 767e7c4fc90SDmitry Karpeev PetscErrorCode DMCreateDecompositionDM(DM dm, const char* name, DM *ddm) 768e7c4fc90SDmitry Karpeev { 769e7c4fc90SDmitry Karpeev PetscErrorCode ierr; 770e7c4fc90SDmitry Karpeev 771e7c4fc90SDmitry Karpeev PetscFunctionBegin; 772e7c4fc90SDmitry Karpeev PetscValidHeaderSpecific(dm,DM_CLASSID,1); 773e7c4fc90SDmitry Karpeev PetscValidCharPointer(name,2); 774e7c4fc90SDmitry Karpeev PetscValidPointer(ddm,3); 775e7c4fc90SDmitry Karpeev if(!dm->ops->createdecompositiondm) { 776e7c4fc90SDmitry Karpeev *ddm = PETSC_NULL; 777e7c4fc90SDmitry Karpeev } 778e7c4fc90SDmitry Karpeev else { 779e7c4fc90SDmitry Karpeev ierr = (*dm->ops->createdecompositiondm)(dm,name,ddm); CHKERRQ(ierr); 780e7c4fc90SDmitry Karpeev } 781e7c4fc90SDmitry Karpeev PetscFunctionReturn(0); 782e7c4fc90SDmitry Karpeev } 783e7c4fc90SDmitry Karpeev 784e7c4fc90SDmitry Karpeev #undef __FUNCT__ 7854d343eeaSMatthew G Knepley #define __FUNCT__ "DMCreateFieldIS" 7864f3b5142SJed Brown /*@C 7874d343eeaSMatthew G Knepley DMCreateFieldIS - Creates a set of IS objects with the global indices of dofs for each field 7884d343eeaSMatthew G Knepley 7894d343eeaSMatthew G Knepley Not collective 7904d343eeaSMatthew G Knepley 7914d343eeaSMatthew G Knepley Input Parameter: 7924d343eeaSMatthew G Knepley . dm - the DM object 7934d343eeaSMatthew G Knepley 7944d343eeaSMatthew G Knepley Output Parameters: 79521c9b008SJed Brown + numFields - The number of fields (or PETSC_NULL if not requested) 79621c9b008SJed Brown . names - The name for each field (or PETSC_NULL if not requested) 79721c9b008SJed Brown - fields - The global indices for each field (or PETSC_NULL if not requested) 7984d343eeaSMatthew G Knepley 7994d343eeaSMatthew G Knepley Level: intermediate 8004d343eeaSMatthew G Knepley 80121c9b008SJed Brown Notes: 80221c9b008SJed Brown The user is responsible for freeing all requested arrays. In particular, every entry of names should be freed with 80321c9b008SJed Brown PetscFree(), every entry of fields should be destroyed with ISDestroy(), and both arrays should be freed with 80421c9b008SJed Brown PetscFree(). 80521c9b008SJed Brown 8064d343eeaSMatthew G Knepley .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix() 8074d343eeaSMatthew G Knepley @*/ 80821c9b008SJed Brown PetscErrorCode DMCreateFieldIS(DM dm, PetscInt *numFields, char ***names, IS **fields) 8094d343eeaSMatthew G Knepley { 8104d343eeaSMatthew G Knepley PetscErrorCode ierr; 8114d343eeaSMatthew G Knepley 8124d343eeaSMatthew G Knepley PetscFunctionBegin; 8134d343eeaSMatthew G Knepley PetscValidHeaderSpecific(dm,DM_CLASSID,1); 8144d343eeaSMatthew G Knepley if (numFields) {PetscValidPointer(numFields,2);} 8154d343eeaSMatthew G Knepley if (names) {PetscValidPointer(names,3);} 8164d343eeaSMatthew G Knepley if (fields) {PetscValidPointer(fields,4);} 8174d343eeaSMatthew G Knepley ierr = (*dm->ops->createfieldis)(dm, numFields, names, fields);CHKERRQ(ierr); 8184d343eeaSMatthew G Knepley PetscFunctionReturn(0); 8194d343eeaSMatthew G Knepley } 8204d343eeaSMatthew G Knepley 8215fe1f584SPeter Brune 822a89ea682SMatthew G Knepley #undef __FUNCT__ 823e7c4fc90SDmitry Karpeev #define __FUNCT__ "DMCreateDecomposition" 824e7c4fc90SDmitry Karpeev /*@C 825e7c4fc90SDmitry Karpeev DMCreateDecomposition - Returns a list of IS objects defining a decomposition of a problem into subproblems: 826e7c4fc90SDmitry Karpeev each IS contains the global indices of the dofs of the corresponding subproblem. 827e7c4fc90SDmitry Karpeev The optional list of DMs define the DM for each subproblem. 828e7c4fc90SDmitry Karpeev Generalizes DMCreateFieldIS(). 829e7c4fc90SDmitry Karpeev 830e7c4fc90SDmitry Karpeev Not collective 831e7c4fc90SDmitry Karpeev 832e7c4fc90SDmitry Karpeev Input Parameter: 833e7c4fc90SDmitry Karpeev . dm - the DM object 834e7c4fc90SDmitry Karpeev 835e7c4fc90SDmitry Karpeev Output Parameters: 836e7c4fc90SDmitry Karpeev + len - The number of subproblems in the decomposition (or PETSC_NULL if not requested) 837e7c4fc90SDmitry Karpeev . namelist - The name for each subproblem (or PETSC_NULL if not requested) 838e7c4fc90SDmitry Karpeev . islist - The global indices for each subproblem (or PETSC_NULL if not requested) 839e7c4fc90SDmitry Karpeev - dmlist - The DMs for each subproblem (or PETSC_NULL, if not requested; if PETSC_NULL is returned, no DMs are defined) 840e7c4fc90SDmitry Karpeev 841e7c4fc90SDmitry Karpeev Level: intermediate 842e7c4fc90SDmitry Karpeev 843e7c4fc90SDmitry Karpeev Notes: 844e7c4fc90SDmitry Karpeev The user is responsible for freeing all requested arrays. In particular, every entry of names should be freed with 845e7c4fc90SDmitry Karpeev PetscFree(), every entry of is should be destroyed with ISDestroy(), every entry of dm should be destroyed with DMDestroy(), 846e7c4fc90SDmitry Karpeev and all of the arrays should be freed with PetscFree(). 847e7c4fc90SDmitry Karpeev 848e7c4fc90SDmitry Karpeev .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMCreateFieldIS() 849e7c4fc90SDmitry Karpeev @*/ 850e7c4fc90SDmitry Karpeev PetscErrorCode DMCreateDecomposition(DM dm, PetscInt *len, char ***namelist, IS **islist, DM **dmlist) 851e7c4fc90SDmitry Karpeev { 852e7c4fc90SDmitry Karpeev PetscErrorCode ierr; 853e7c4fc90SDmitry Karpeev 854e7c4fc90SDmitry Karpeev PetscFunctionBegin; 855e7c4fc90SDmitry Karpeev PetscValidHeaderSpecific(dm,DM_CLASSID,1); 856e7c4fc90SDmitry Karpeev if (len) {PetscValidPointer(len,2);} 857e7c4fc90SDmitry Karpeev if (namelist) {PetscValidPointer(namelist,3);} 858e7c4fc90SDmitry Karpeev if (islist) {PetscValidPointer(islist,4);} 859e7c4fc90SDmitry Karpeev if (dmlist) {PetscValidPointer(dmlist,5);} 860e7c4fc90SDmitry Karpeev if(!dm->ops->createdecomposition) { 861e7c4fc90SDmitry Karpeev ierr = (*dm->ops->createfieldis)(dm, len, namelist, islist);CHKERRQ(ierr); 862e7c4fc90SDmitry Karpeev /* By default there are no DMs associated with subproblems. */ 863e7c4fc90SDmitry Karpeev if(dmlist) *dmlist = PETSC_NULL; 864e7c4fc90SDmitry Karpeev } 865e7c4fc90SDmitry Karpeev else { 866e7c4fc90SDmitry Karpeev ierr = (*dm->ops->createdecomposition)(dm,len,namelist,islist,dmlist); CHKERRQ(ierr); 867e7c4fc90SDmitry Karpeev } 868e7c4fc90SDmitry Karpeev PetscFunctionReturn(0); 869e7c4fc90SDmitry Karpeev } 870e7c4fc90SDmitry Karpeev 871e7c4fc90SDmitry Karpeev 872e7c4fc90SDmitry Karpeev #undef __FUNCT__ 87347c6ae99SBarry Smith #define __FUNCT__ "DMRefine" 87447c6ae99SBarry Smith /*@ 87547c6ae99SBarry Smith DMRefine - Refines a DM object 87647c6ae99SBarry Smith 87747c6ae99SBarry Smith Collective on DM 87847c6ae99SBarry Smith 87947c6ae99SBarry Smith Input Parameter: 88047c6ae99SBarry Smith + dm - the DM object 88191d95f02SJed Brown - comm - the communicator to contain the new DM object (or MPI_COMM_NULL) 88247c6ae99SBarry Smith 88347c6ae99SBarry Smith Output Parameter: 884ae0a1c52SMatthew G Knepley . dmf - the refined DM, or PETSC_NULL 885ae0a1c52SMatthew G Knepley 886ae0a1c52SMatthew G Knepley Note: If no refinement was done, the return value is PETSC_NULL 88747c6ae99SBarry Smith 88847c6ae99SBarry Smith Level: developer 88947c6ae99SBarry Smith 890e727c939SJed Brown .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation() 89147c6ae99SBarry Smith @*/ 8927087cfbeSBarry Smith PetscErrorCode DMRefine(DM dm,MPI_Comm comm,DM *dmf) 89347c6ae99SBarry Smith { 89447c6ae99SBarry Smith PetscErrorCode ierr; 89547c6ae99SBarry Smith 89647c6ae99SBarry Smith PetscFunctionBegin; 897732e2eb9SMatthew G Knepley PetscValidHeaderSpecific(dm,DM_CLASSID,1); 89847c6ae99SBarry Smith ierr = (*dm->ops->refine)(dm,comm,dmf);CHKERRQ(ierr); 8994057135bSMatthew G Knepley if (*dmf) { 90043842a1eSJed Brown (*dmf)->ops->creatematrix = dm->ops->creatematrix; 901644e2e5bSBarry Smith (*dmf)->ops->initialguess = dm->ops->initialguess; 902644e2e5bSBarry Smith (*dmf)->ops->function = dm->ops->function; 903644e2e5bSBarry Smith (*dmf)->ops->functionj = dm->ops->functionj; 904644e2e5bSBarry Smith if (dm->ops->jacobian != DMComputeJacobianDefault) { 905644e2e5bSBarry Smith (*dmf)->ops->jacobian = dm->ops->jacobian; 906644e2e5bSBarry Smith } 9078cd211a4SJed Brown ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)dm,(PetscObject)*dmf);CHKERRQ(ierr); 908644e2e5bSBarry Smith (*dmf)->ctx = dm->ctx; 909656b349aSBarry Smith (*dmf)->levelup = dm->levelup + 1; 9104057135bSMatthew G Knepley } 91147c6ae99SBarry Smith PetscFunctionReturn(0); 91247c6ae99SBarry Smith } 91347c6ae99SBarry Smith 91447c6ae99SBarry Smith #undef __FUNCT__ 915eb3f98d2SBarry Smith #define __FUNCT__ "DMGetRefineLevel" 916eb3f98d2SBarry Smith /*@ 917eb3f98d2SBarry Smith DMGetRefineLevel - Get's the number of refinements that have generated this DM. 918eb3f98d2SBarry Smith 919eb3f98d2SBarry Smith Not Collective 920eb3f98d2SBarry Smith 921eb3f98d2SBarry Smith Input Parameter: 922eb3f98d2SBarry Smith . dm - the DM object 923eb3f98d2SBarry Smith 924eb3f98d2SBarry Smith Output Parameter: 925eb3f98d2SBarry Smith . level - number of refinements 926eb3f98d2SBarry Smith 927eb3f98d2SBarry Smith Level: developer 928eb3f98d2SBarry Smith 9296a7d9d85SPeter Brune .seealso DMCoarsen(), DMGetCoarsenLevel(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation() 930eb3f98d2SBarry Smith 931eb3f98d2SBarry Smith @*/ 932eb3f98d2SBarry Smith PetscErrorCode DMGetRefineLevel(DM dm,PetscInt *level) 933eb3f98d2SBarry Smith { 934eb3f98d2SBarry Smith PetscFunctionBegin; 935eb3f98d2SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 936eb3f98d2SBarry Smith *level = dm->levelup; 937eb3f98d2SBarry Smith PetscFunctionReturn(0); 938eb3f98d2SBarry Smith } 939eb3f98d2SBarry Smith 940eb3f98d2SBarry Smith #undef __FUNCT__ 94147c6ae99SBarry Smith #define __FUNCT__ "DMGlobalToLocalBegin" 94247c6ae99SBarry Smith /*@ 94347c6ae99SBarry Smith DMGlobalToLocalBegin - Begins updating local vectors from global vector 94447c6ae99SBarry Smith 94547c6ae99SBarry Smith Neighbor-wise Collective on DM 94647c6ae99SBarry Smith 94747c6ae99SBarry Smith Input Parameters: 94847c6ae99SBarry Smith + dm - the DM object 94947c6ae99SBarry Smith . g - the global vector 95047c6ae99SBarry Smith . mode - INSERT_VALUES or ADD_VALUES 95147c6ae99SBarry Smith - l - the local vector 95247c6ae99SBarry Smith 95347c6ae99SBarry Smith 95447c6ae99SBarry Smith Level: beginner 95547c6ae99SBarry Smith 956e727c939SJed Brown .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin() 95747c6ae99SBarry Smith 95847c6ae99SBarry Smith @*/ 9597087cfbeSBarry Smith PetscErrorCode DMGlobalToLocalBegin(DM dm,Vec g,InsertMode mode,Vec l) 96047c6ae99SBarry Smith { 96147c6ae99SBarry Smith PetscErrorCode ierr; 96247c6ae99SBarry Smith 96347c6ae99SBarry Smith PetscFunctionBegin; 964171400e9SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 965843c4018SMatthew G Knepley ierr = (*dm->ops->globaltolocalbegin)(dm,g,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),l);CHKERRQ(ierr); 96647c6ae99SBarry Smith PetscFunctionReturn(0); 96747c6ae99SBarry Smith } 96847c6ae99SBarry Smith 96947c6ae99SBarry Smith #undef __FUNCT__ 97047c6ae99SBarry Smith #define __FUNCT__ "DMGlobalToLocalEnd" 97147c6ae99SBarry Smith /*@ 97247c6ae99SBarry Smith DMGlobalToLocalEnd - Ends updating local vectors from global vector 97347c6ae99SBarry Smith 97447c6ae99SBarry Smith Neighbor-wise Collective on DM 97547c6ae99SBarry Smith 97647c6ae99SBarry Smith Input Parameters: 97747c6ae99SBarry Smith + dm - the DM object 97847c6ae99SBarry Smith . g - the global vector 97947c6ae99SBarry Smith . mode - INSERT_VALUES or ADD_VALUES 98047c6ae99SBarry Smith - l - the local vector 98147c6ae99SBarry Smith 98247c6ae99SBarry Smith 98347c6ae99SBarry Smith Level: beginner 98447c6ae99SBarry Smith 985e727c939SJed Brown .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin() 98647c6ae99SBarry Smith 98747c6ae99SBarry Smith @*/ 9887087cfbeSBarry Smith PetscErrorCode DMGlobalToLocalEnd(DM dm,Vec g,InsertMode mode,Vec l) 98947c6ae99SBarry Smith { 99047c6ae99SBarry Smith PetscErrorCode ierr; 99147c6ae99SBarry Smith 99247c6ae99SBarry Smith PetscFunctionBegin; 993171400e9SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 994843c4018SMatthew G Knepley ierr = (*dm->ops->globaltolocalend)(dm,g,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),l);CHKERRQ(ierr); 99547c6ae99SBarry Smith PetscFunctionReturn(0); 99647c6ae99SBarry Smith } 99747c6ae99SBarry Smith 99847c6ae99SBarry Smith #undef __FUNCT__ 9999a42bb27SBarry Smith #define __FUNCT__ "DMLocalToGlobalBegin" 100047c6ae99SBarry Smith /*@ 10019a42bb27SBarry Smith DMLocalToGlobalBegin - updates global vectors from local vectors 10029a42bb27SBarry Smith 10039a42bb27SBarry Smith Neighbor-wise Collective on DM 10049a42bb27SBarry Smith 10059a42bb27SBarry Smith Input Parameters: 10069a42bb27SBarry Smith + dm - the DM object 1007f6813fd5SJed Brown . l - the local vector 10089a42bb27SBarry Smith . mode - if INSERT_VALUES then no parallel communication is used, if ADD_VALUES then all ghost points from the same base point accumulate into that 10099a42bb27SBarry Smith base point. 1010f6813fd5SJed Brown - - the global vector 10119a42bb27SBarry Smith 10129a42bb27SBarry Smith Notes: In the ADD_VALUES case you normally would zero the receiving vector before beginning this operation. If you would like to simply add the non-ghosted values in the local 10139a42bb27SBarry Smith array into the global array you need to either (1) zero the ghosted locations and use ADD_VALUES or (2) use INSERT_VALUES into a work global array and then add the work 10149a42bb27SBarry Smith global array to the final global array with VecAXPY(). 10159a42bb27SBarry Smith 10169a42bb27SBarry Smith Level: beginner 10179a42bb27SBarry Smith 1018e727c939SJed Brown .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalBegin() 10199a42bb27SBarry Smith 10209a42bb27SBarry Smith @*/ 10217087cfbeSBarry Smith PetscErrorCode DMLocalToGlobalBegin(DM dm,Vec l,InsertMode mode,Vec g) 10229a42bb27SBarry Smith { 10239a42bb27SBarry Smith PetscErrorCode ierr; 10249a42bb27SBarry Smith 10259a42bb27SBarry Smith PetscFunctionBegin; 1026171400e9SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1027843c4018SMatthew G Knepley ierr = (*dm->ops->localtoglobalbegin)(dm,l,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),g);CHKERRQ(ierr); 10289a42bb27SBarry Smith PetscFunctionReturn(0); 10299a42bb27SBarry Smith } 10309a42bb27SBarry Smith 10319a42bb27SBarry Smith #undef __FUNCT__ 10329a42bb27SBarry Smith #define __FUNCT__ "DMLocalToGlobalEnd" 10339a42bb27SBarry Smith /*@ 10349a42bb27SBarry Smith DMLocalToGlobalEnd - updates global vectors from local vectors 103547c6ae99SBarry Smith 103647c6ae99SBarry Smith Neighbor-wise Collective on DM 103747c6ae99SBarry Smith 103847c6ae99SBarry Smith Input Parameters: 103947c6ae99SBarry Smith + dm - the DM object 1040f6813fd5SJed Brown . l - the local vector 104147c6ae99SBarry Smith . mode - INSERT_VALUES or ADD_VALUES 1042f6813fd5SJed Brown - g - the global vector 104347c6ae99SBarry Smith 104447c6ae99SBarry Smith 104547c6ae99SBarry Smith Level: beginner 104647c6ae99SBarry Smith 1047e727c939SJed Brown .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalEnd() 104847c6ae99SBarry Smith 104947c6ae99SBarry Smith @*/ 10507087cfbeSBarry Smith PetscErrorCode DMLocalToGlobalEnd(DM dm,Vec l,InsertMode mode,Vec g) 105147c6ae99SBarry Smith { 105247c6ae99SBarry Smith PetscErrorCode ierr; 105347c6ae99SBarry Smith 105447c6ae99SBarry Smith PetscFunctionBegin; 1055171400e9SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1056843c4018SMatthew G Knepley ierr = (*dm->ops->localtoglobalend)(dm,l,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),g);CHKERRQ(ierr); 105747c6ae99SBarry Smith PetscFunctionReturn(0); 105847c6ae99SBarry Smith } 105947c6ae99SBarry Smith 106047c6ae99SBarry Smith #undef __FUNCT__ 106147c6ae99SBarry Smith #define __FUNCT__ "DMComputeJacobianDefault" 106247c6ae99SBarry Smith /*@ 106347c6ae99SBarry Smith DMComputeJacobianDefault - computes the Jacobian using the DMComputeFunction() if Jacobian computer is not provided 106447c6ae99SBarry Smith 106547c6ae99SBarry Smith Collective on DM 106647c6ae99SBarry Smith 106747c6ae99SBarry Smith Input Parameter: 106847c6ae99SBarry Smith + dm - the DM object 106947c6ae99SBarry Smith . x - location to compute Jacobian at; may be ignored for linear problems 107047c6ae99SBarry Smith . A - matrix that defines the operator for the linear solve 107147c6ae99SBarry Smith - B - the matrix used to construct the preconditioner 107247c6ae99SBarry Smith 107347c6ae99SBarry Smith Level: developer 107447c6ae99SBarry Smith 1075e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 107647c6ae99SBarry Smith DMSetFunction() 107747c6ae99SBarry Smith 107847c6ae99SBarry Smith @*/ 10797087cfbeSBarry Smith PetscErrorCode DMComputeJacobianDefault(DM dm,Vec x,Mat A,Mat B,MatStructure *stflag) 108047c6ae99SBarry Smith { 108147c6ae99SBarry Smith PetscErrorCode ierr; 1082171400e9SBarry Smith 108347c6ae99SBarry Smith PetscFunctionBegin; 1084171400e9SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 108547c6ae99SBarry Smith *stflag = SAME_NONZERO_PATTERN; 108647c6ae99SBarry Smith ierr = MatFDColoringApply(B,dm->fd,x,stflag,dm);CHKERRQ(ierr); 108747c6ae99SBarry Smith if (A != B) { 108847c6ae99SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 108947c6ae99SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 109047c6ae99SBarry Smith } 109147c6ae99SBarry Smith PetscFunctionReturn(0); 109247c6ae99SBarry Smith } 109347c6ae99SBarry Smith 109447c6ae99SBarry Smith #undef __FUNCT__ 109547c6ae99SBarry Smith #define __FUNCT__ "DMCoarsen" 109647c6ae99SBarry Smith /*@ 109747c6ae99SBarry Smith DMCoarsen - Coarsens a DM object 109847c6ae99SBarry Smith 109947c6ae99SBarry Smith Collective on DM 110047c6ae99SBarry Smith 110147c6ae99SBarry Smith Input Parameter: 110247c6ae99SBarry Smith + dm - the DM object 110391d95f02SJed Brown - comm - the communicator to contain the new DM object (or MPI_COMM_NULL) 110447c6ae99SBarry Smith 110547c6ae99SBarry Smith Output Parameter: 110647c6ae99SBarry Smith . dmc - the coarsened DM 110747c6ae99SBarry Smith 110847c6ae99SBarry Smith Level: developer 110947c6ae99SBarry Smith 1110e727c939SJed Brown .seealso DMRefine(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation() 111147c6ae99SBarry Smith 111247c6ae99SBarry Smith @*/ 11137087cfbeSBarry Smith PetscErrorCode DMCoarsen(DM dm, MPI_Comm comm, DM *dmc) 111447c6ae99SBarry Smith { 111547c6ae99SBarry Smith PetscErrorCode ierr; 1116b17ce1afSJed Brown DMCoarsenHookLink link; 111747c6ae99SBarry Smith 111847c6ae99SBarry Smith PetscFunctionBegin; 1119171400e9SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 112047c6ae99SBarry Smith ierr = (*dm->ops->coarsen)(dm, comm, dmc);CHKERRQ(ierr); 112143842a1eSJed Brown (*dmc)->ops->creatematrix = dm->ops->creatematrix; 112247c6ae99SBarry Smith (*dmc)->ops->initialguess = dm->ops->initialguess; 112347c6ae99SBarry Smith (*dmc)->ops->function = dm->ops->function; 112447c6ae99SBarry Smith (*dmc)->ops->functionj = dm->ops->functionj; 112547c6ae99SBarry Smith if (dm->ops->jacobian != DMComputeJacobianDefault) { 112647c6ae99SBarry Smith (*dmc)->ops->jacobian = dm->ops->jacobian; 112747c6ae99SBarry Smith } 11288cd211a4SJed Brown ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)dm,(PetscObject)*dmc);CHKERRQ(ierr); 1129644e2e5bSBarry Smith (*dmc)->ctx = dm->ctx; 1130656b349aSBarry Smith (*dmc)->leveldown = dm->leveldown + 1; 1131b17ce1afSJed Brown for (link=dm->coarsenhook; link; link=link->next) { 1132b17ce1afSJed Brown if (link->coarsenhook) {ierr = (*link->coarsenhook)(dm,*dmc,link->ctx);CHKERRQ(ierr);} 1133b17ce1afSJed Brown } 1134b17ce1afSJed Brown PetscFunctionReturn(0); 1135b17ce1afSJed Brown } 1136b17ce1afSJed Brown 1137b17ce1afSJed Brown #undef __FUNCT__ 1138b17ce1afSJed Brown #define __FUNCT__ "DMCoarsenHookAdd" 1139b17ce1afSJed Brown /*@ 1140b17ce1afSJed Brown DMCoarsenHookAdd - adds a callback to be run when restricting a nonlinear problem to the coarse grid 1141b17ce1afSJed Brown 1142b17ce1afSJed Brown Logically Collective 1143b17ce1afSJed Brown 1144b17ce1afSJed Brown Input Arguments: 1145b17ce1afSJed Brown + fine - nonlinear solver context on which to run a hook when restricting to a coarser level 1146b17ce1afSJed Brown . coarsenhook - function to run when setting up a coarser level 1147b17ce1afSJed Brown . restricthook - function to run to update data on coarser levels (once per SNESSolve()) 1148b17ce1afSJed Brown - ctx - [optional] user-defined context for provide data for the hooks (may be PETSC_NULL) 1149b17ce1afSJed Brown 1150b17ce1afSJed Brown Calling sequence of coarsenhook: 1151b17ce1afSJed Brown $ coarsenhook(DM fine,DM coarse,void *ctx); 1152b17ce1afSJed Brown 1153b17ce1afSJed Brown + fine - fine level DM 1154b17ce1afSJed Brown . coarse - coarse level DM to restrict problem to 1155b17ce1afSJed Brown - ctx - optional user-defined function context 1156b17ce1afSJed Brown 1157b17ce1afSJed Brown Calling sequence for restricthook: 1158b17ce1afSJed Brown $ restricthook(DM fine,Mat mrestrict,Mat inject,DM coarse,void *ctx) 1159b17ce1afSJed Brown 1160b17ce1afSJed Brown + fine - fine level DM 1161b17ce1afSJed Brown . mrestrict - matrix restricting a fine-level solution to the coarse grid 1162b17ce1afSJed Brown . inject - matrix restricting by applying the transpose of injection 1163b17ce1afSJed Brown . coarse - coarse level DM to update 1164b17ce1afSJed Brown - ctx - optional user-defined function context 1165b17ce1afSJed Brown 1166b17ce1afSJed Brown Level: advanced 1167b17ce1afSJed Brown 1168b17ce1afSJed Brown Notes: 1169b17ce1afSJed Brown This function is only needed if auxiliary data needs to be set up on coarse grids. 1170b17ce1afSJed Brown 1171b17ce1afSJed Brown If this function is called multiple times, the hooks will be run in the order they are added. 1172b17ce1afSJed Brown 1173b17ce1afSJed Brown The 1174b17ce1afSJed Brown 1175b17ce1afSJed Brown In order to compose with nonlinear preconditioning without duplicating storage, the hook should be implemented to 1176b17ce1afSJed Brown extract the finest level information from its context (instead of from the SNES). 1177b17ce1afSJed Brown 1178b17ce1afSJed Brown .seealso: SNESFASGetInterpolation(), SNESFASGetInjection(), PetscObjectCompose(), PetscContainerCreate() 1179b17ce1afSJed Brown @*/ 1180b17ce1afSJed Brown PetscErrorCode DMCoarsenHookAdd(DM fine,PetscErrorCode (*coarsenhook)(DM,DM,void*),PetscErrorCode (*restricthook)(DM,Mat,Vec,Mat,DM,void*),void *ctx) 1181b17ce1afSJed Brown { 1182b17ce1afSJed Brown PetscErrorCode ierr; 1183b17ce1afSJed Brown DMCoarsenHookLink link,*p; 1184b17ce1afSJed Brown 1185b17ce1afSJed Brown PetscFunctionBegin; 1186b17ce1afSJed Brown PetscValidHeaderSpecific(fine,DM_CLASSID,1); 11876bfea28cSJed Brown for (p=&fine->coarsenhook; *p; p=&(*p)->next) {} /* Scan to the end of the current list of hooks */ 1188b17ce1afSJed Brown ierr = PetscMalloc(sizeof(struct _DMCoarsenHookLink),&link);CHKERRQ(ierr); 1189b17ce1afSJed Brown link->coarsenhook = coarsenhook; 1190b17ce1afSJed Brown link->restricthook = restricthook; 1191b17ce1afSJed Brown link->ctx = ctx; 11926cab3a1bSJed Brown link->next = PETSC_NULL; 1193b17ce1afSJed Brown *p = link; 1194b17ce1afSJed Brown PetscFunctionReturn(0); 1195b17ce1afSJed Brown } 1196b17ce1afSJed Brown 1197b17ce1afSJed Brown #undef __FUNCT__ 1198b17ce1afSJed Brown #define __FUNCT__ "DMRestrict" 1199b17ce1afSJed Brown /*@ 1200b17ce1afSJed Brown DMRestrict - restricts user-defined problem data to a coarser DM by running hooks registered by DMCoarsenHookAdd() 1201b17ce1afSJed Brown 1202b17ce1afSJed Brown Collective if any hooks are 1203b17ce1afSJed Brown 1204b17ce1afSJed Brown Input Arguments: 1205b17ce1afSJed Brown + fine - finer DM to use as a base 1206b17ce1afSJed Brown . restrct - restriction matrix, apply using MatRestrict() 1207b17ce1afSJed Brown . inject - injection matrix, also use MatRestrict() 1208b17ce1afSJed Brown - coarse - coarer DM to update 1209b17ce1afSJed Brown 1210b17ce1afSJed Brown Level: developer 1211b17ce1afSJed Brown 1212b17ce1afSJed Brown .seealso: DMCoarsenHookAdd(), MatRestrict() 1213b17ce1afSJed Brown @*/ 1214b17ce1afSJed Brown PetscErrorCode DMRestrict(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse) 1215b17ce1afSJed Brown { 1216b17ce1afSJed Brown PetscErrorCode ierr; 1217b17ce1afSJed Brown DMCoarsenHookLink link; 1218b17ce1afSJed Brown 1219b17ce1afSJed Brown PetscFunctionBegin; 1220b17ce1afSJed Brown for (link=fine->coarsenhook; link; link=link->next) { 1221b17ce1afSJed Brown if (link->restricthook) {ierr = (*link->restricthook)(fine,restrct,rscale,inject,coarse,link->ctx);CHKERRQ(ierr);} 1222b17ce1afSJed Brown } 122347c6ae99SBarry Smith PetscFunctionReturn(0); 122447c6ae99SBarry Smith } 122547c6ae99SBarry Smith 122647c6ae99SBarry Smith #undef __FUNCT__ 12275fe1f584SPeter Brune #define __FUNCT__ "DMGetCoarsenLevel" 12285fe1f584SPeter Brune /*@ 12296a7d9d85SPeter Brune DMGetCoarsenLevel - Get's the number of coarsenings that have generated this DM. 12305fe1f584SPeter Brune 12315fe1f584SPeter Brune Not Collective 12325fe1f584SPeter Brune 12335fe1f584SPeter Brune Input Parameter: 12345fe1f584SPeter Brune . dm - the DM object 12355fe1f584SPeter Brune 12365fe1f584SPeter Brune Output Parameter: 12376a7d9d85SPeter Brune . level - number of coarsenings 12385fe1f584SPeter Brune 12395fe1f584SPeter Brune Level: developer 12405fe1f584SPeter Brune 12416a7d9d85SPeter Brune .seealso DMCoarsen(), DMGetRefineLevel(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation() 12425fe1f584SPeter Brune 12435fe1f584SPeter Brune @*/ 12445fe1f584SPeter Brune PetscErrorCode DMGetCoarsenLevel(DM dm,PetscInt *level) 12455fe1f584SPeter Brune { 12465fe1f584SPeter Brune PetscFunctionBegin; 12475fe1f584SPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 12485fe1f584SPeter Brune *level = dm->leveldown; 12495fe1f584SPeter Brune PetscFunctionReturn(0); 12505fe1f584SPeter Brune } 12515fe1f584SPeter Brune 12525fe1f584SPeter Brune 12535fe1f584SPeter Brune 12545fe1f584SPeter Brune #undef __FUNCT__ 125547c6ae99SBarry Smith #define __FUNCT__ "DMRefineHierarchy" 125647c6ae99SBarry Smith /*@C 125747c6ae99SBarry Smith DMRefineHierarchy - Refines a DM object, all levels at once 125847c6ae99SBarry Smith 125947c6ae99SBarry Smith Collective on DM 126047c6ae99SBarry Smith 126147c6ae99SBarry Smith Input Parameter: 126247c6ae99SBarry Smith + dm - the DM object 126347c6ae99SBarry Smith - nlevels - the number of levels of refinement 126447c6ae99SBarry Smith 126547c6ae99SBarry Smith Output Parameter: 126647c6ae99SBarry Smith . dmf - the refined DM hierarchy 126747c6ae99SBarry Smith 126847c6ae99SBarry Smith Level: developer 126947c6ae99SBarry Smith 1270e727c939SJed Brown .seealso DMCoarsenHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation() 127147c6ae99SBarry Smith 127247c6ae99SBarry Smith @*/ 12737087cfbeSBarry Smith PetscErrorCode DMRefineHierarchy(DM dm,PetscInt nlevels,DM dmf[]) 127447c6ae99SBarry Smith { 127547c6ae99SBarry Smith PetscErrorCode ierr; 127647c6ae99SBarry Smith 127747c6ae99SBarry Smith PetscFunctionBegin; 1278171400e9SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 127947c6ae99SBarry Smith if (nlevels < 0) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative"); 128047c6ae99SBarry Smith if (nlevels == 0) PetscFunctionReturn(0); 128147c6ae99SBarry Smith if (dm->ops->refinehierarchy) { 128247c6ae99SBarry Smith ierr = (*dm->ops->refinehierarchy)(dm,nlevels,dmf);CHKERRQ(ierr); 128347c6ae99SBarry Smith } else if (dm->ops->refine) { 128447c6ae99SBarry Smith PetscInt i; 128547c6ae99SBarry Smith 128647c6ae99SBarry Smith ierr = DMRefine(dm,((PetscObject)dm)->comm,&dmf[0]);CHKERRQ(ierr); 128747c6ae99SBarry Smith for (i=1; i<nlevels; i++) { 128847c6ae99SBarry Smith ierr = DMRefine(dmf[i-1],((PetscObject)dm)->comm,&dmf[i]);CHKERRQ(ierr); 128947c6ae99SBarry Smith } 129047c6ae99SBarry Smith } else { 129147c6ae99SBarry Smith SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No RefineHierarchy for this DM yet"); 129247c6ae99SBarry Smith } 129347c6ae99SBarry Smith PetscFunctionReturn(0); 129447c6ae99SBarry Smith } 129547c6ae99SBarry Smith 129647c6ae99SBarry Smith #undef __FUNCT__ 129747c6ae99SBarry Smith #define __FUNCT__ "DMCoarsenHierarchy" 129847c6ae99SBarry Smith /*@C 129947c6ae99SBarry Smith DMCoarsenHierarchy - Coarsens a DM object, all levels at once 130047c6ae99SBarry Smith 130147c6ae99SBarry Smith Collective on DM 130247c6ae99SBarry Smith 130347c6ae99SBarry Smith Input Parameter: 130447c6ae99SBarry Smith + dm - the DM object 130547c6ae99SBarry Smith - nlevels - the number of levels of coarsening 130647c6ae99SBarry Smith 130747c6ae99SBarry Smith Output Parameter: 130847c6ae99SBarry Smith . dmc - the coarsened DM hierarchy 130947c6ae99SBarry Smith 131047c6ae99SBarry Smith Level: developer 131147c6ae99SBarry Smith 1312e727c939SJed Brown .seealso DMRefineHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation() 131347c6ae99SBarry Smith 131447c6ae99SBarry Smith @*/ 13157087cfbeSBarry Smith PetscErrorCode DMCoarsenHierarchy(DM dm, PetscInt nlevels, DM dmc[]) 131647c6ae99SBarry Smith { 131747c6ae99SBarry Smith PetscErrorCode ierr; 131847c6ae99SBarry Smith 131947c6ae99SBarry Smith PetscFunctionBegin; 1320171400e9SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 132147c6ae99SBarry Smith if (nlevels < 0) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative"); 132247c6ae99SBarry Smith if (nlevels == 0) PetscFunctionReturn(0); 132347c6ae99SBarry Smith PetscValidPointer(dmc,3); 132447c6ae99SBarry Smith if (dm->ops->coarsenhierarchy) { 132547c6ae99SBarry Smith ierr = (*dm->ops->coarsenhierarchy)(dm, nlevels, dmc);CHKERRQ(ierr); 132647c6ae99SBarry Smith } else if (dm->ops->coarsen) { 132747c6ae99SBarry Smith PetscInt i; 132847c6ae99SBarry Smith 132947c6ae99SBarry Smith ierr = DMCoarsen(dm,((PetscObject)dm)->comm,&dmc[0]);CHKERRQ(ierr); 133047c6ae99SBarry Smith for (i=1; i<nlevels; i++) { 133147c6ae99SBarry Smith ierr = DMCoarsen(dmc[i-1],((PetscObject)dm)->comm,&dmc[i]);CHKERRQ(ierr); 133247c6ae99SBarry Smith } 133347c6ae99SBarry Smith } else { 133447c6ae99SBarry Smith SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No CoarsenHierarchy for this DM yet"); 133547c6ae99SBarry Smith } 133647c6ae99SBarry Smith PetscFunctionReturn(0); 133747c6ae99SBarry Smith } 133847c6ae99SBarry Smith 133947c6ae99SBarry Smith #undef __FUNCT__ 1340e727c939SJed Brown #define __FUNCT__ "DMCreateAggregates" 134147c6ae99SBarry Smith /*@ 1342e727c939SJed Brown DMCreateAggregates - Gets the aggregates that map between 134347c6ae99SBarry Smith grids associated with two DMs. 134447c6ae99SBarry Smith 134547c6ae99SBarry Smith Collective on DM 134647c6ae99SBarry Smith 134747c6ae99SBarry Smith Input Parameters: 134847c6ae99SBarry Smith + dmc - the coarse grid DM 134947c6ae99SBarry Smith - dmf - the fine grid DM 135047c6ae99SBarry Smith 135147c6ae99SBarry Smith Output Parameters: 135247c6ae99SBarry Smith . rest - the restriction matrix (transpose of the projection matrix) 135347c6ae99SBarry Smith 135447c6ae99SBarry Smith Level: intermediate 135547c6ae99SBarry Smith 135647c6ae99SBarry Smith .keywords: interpolation, restriction, multigrid 135747c6ae99SBarry Smith 1358e727c939SJed Brown .seealso: DMRefine(), DMCreateInjection(), DMCreateInterpolation() 135947c6ae99SBarry Smith @*/ 1360e727c939SJed Brown PetscErrorCode DMCreateAggregates(DM dmc, DM dmf, Mat *rest) 136147c6ae99SBarry Smith { 136247c6ae99SBarry Smith PetscErrorCode ierr; 136347c6ae99SBarry Smith 136447c6ae99SBarry Smith PetscFunctionBegin; 1365171400e9SBarry Smith PetscValidHeaderSpecific(dmc,DM_CLASSID,1); 1366171400e9SBarry Smith PetscValidHeaderSpecific(dmf,DM_CLASSID,2); 136747c6ae99SBarry Smith ierr = (*dmc->ops->getaggregates)(dmc, dmf, rest);CHKERRQ(ierr); 136847c6ae99SBarry Smith PetscFunctionReturn(0); 136947c6ae99SBarry Smith } 137047c6ae99SBarry Smith 137147c6ae99SBarry Smith #undef __FUNCT__ 13721a266240SBarry Smith #define __FUNCT__ "DMSetApplicationContextDestroy" 13731a266240SBarry Smith /*@C 13741a266240SBarry Smith DMSetApplicationContextDestroy - Sets a user function that will be called to destroy the application context when the DM is destroyed 13751a266240SBarry Smith 13761a266240SBarry Smith Not Collective 13771a266240SBarry Smith 13781a266240SBarry Smith Input Parameters: 13791a266240SBarry Smith + dm - the DM object 13801a266240SBarry Smith - destroy - the destroy function 13811a266240SBarry Smith 13821a266240SBarry Smith Level: intermediate 13831a266240SBarry Smith 1384e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext() 13851a266240SBarry Smith 1386f07f9ceaSJed Brown @*/ 13871a266240SBarry Smith PetscErrorCode DMSetApplicationContextDestroy(DM dm,PetscErrorCode (*destroy)(void**)) 13881a266240SBarry Smith { 13891a266240SBarry Smith PetscFunctionBegin; 1390171400e9SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 13911a266240SBarry Smith dm->ctxdestroy = destroy; 13921a266240SBarry Smith PetscFunctionReturn(0); 13931a266240SBarry Smith } 13941a266240SBarry Smith 13951a266240SBarry Smith #undef __FUNCT__ 13961b2093e4SBarry Smith #define __FUNCT__ "DMSetApplicationContext" 1397b07ff414SBarry Smith /*@ 13981b2093e4SBarry Smith DMSetApplicationContext - Set a user context into a DM object 139947c6ae99SBarry Smith 140047c6ae99SBarry Smith Not Collective 140147c6ae99SBarry Smith 140247c6ae99SBarry Smith Input Parameters: 140347c6ae99SBarry Smith + dm - the DM object 140447c6ae99SBarry Smith - ctx - the user context 140547c6ae99SBarry Smith 140647c6ae99SBarry Smith Level: intermediate 140747c6ae99SBarry Smith 1408e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext() 140947c6ae99SBarry Smith 141047c6ae99SBarry Smith @*/ 14111b2093e4SBarry Smith PetscErrorCode DMSetApplicationContext(DM dm,void *ctx) 141247c6ae99SBarry Smith { 141347c6ae99SBarry Smith PetscFunctionBegin; 1414171400e9SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 141547c6ae99SBarry Smith dm->ctx = ctx; 141647c6ae99SBarry Smith PetscFunctionReturn(0); 141747c6ae99SBarry Smith } 141847c6ae99SBarry Smith 141947c6ae99SBarry Smith #undef __FUNCT__ 14201b2093e4SBarry Smith #define __FUNCT__ "DMGetApplicationContext" 142147c6ae99SBarry Smith /*@ 14221b2093e4SBarry Smith DMGetApplicationContext - Gets a user context from a DM object 142347c6ae99SBarry Smith 142447c6ae99SBarry Smith Not Collective 142547c6ae99SBarry Smith 142647c6ae99SBarry Smith Input Parameter: 142747c6ae99SBarry Smith . dm - the DM object 142847c6ae99SBarry Smith 142947c6ae99SBarry Smith Output Parameter: 143047c6ae99SBarry Smith . ctx - the user context 143147c6ae99SBarry Smith 143247c6ae99SBarry Smith Level: intermediate 143347c6ae99SBarry Smith 1434e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext() 143547c6ae99SBarry Smith 143647c6ae99SBarry Smith @*/ 14371b2093e4SBarry Smith PetscErrorCode DMGetApplicationContext(DM dm,void *ctx) 143847c6ae99SBarry Smith { 143947c6ae99SBarry Smith PetscFunctionBegin; 1440171400e9SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 14411b2093e4SBarry Smith *(void**)ctx = dm->ctx; 144247c6ae99SBarry Smith PetscFunctionReturn(0); 144347c6ae99SBarry Smith } 144447c6ae99SBarry Smith 144547c6ae99SBarry Smith #undef __FUNCT__ 144647c6ae99SBarry Smith #define __FUNCT__ "DMSetInitialGuess" 14477e833e3aSBarry Smith /*@C 144847c6ae99SBarry Smith DMSetInitialGuess - sets a function to compute an initial guess vector entries for the solvers 144947c6ae99SBarry Smith 145047c6ae99SBarry Smith Logically Collective on DM 145147c6ae99SBarry Smith 145247c6ae99SBarry Smith Input Parameter: 145347c6ae99SBarry Smith + dm - the DM object to destroy 145447c6ae99SBarry Smith - f - the function to compute the initial guess 145547c6ae99SBarry Smith 145647c6ae99SBarry Smith Level: intermediate 145747c6ae99SBarry Smith 1458e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian() 145947c6ae99SBarry Smith 1460f07f9ceaSJed Brown @*/ 14617087cfbeSBarry Smith PetscErrorCode DMSetInitialGuess(DM dm,PetscErrorCode (*f)(DM,Vec)) 146247c6ae99SBarry Smith { 146347c6ae99SBarry Smith PetscFunctionBegin; 1464171400e9SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 146547c6ae99SBarry Smith dm->ops->initialguess = f; 146647c6ae99SBarry Smith PetscFunctionReturn(0); 146747c6ae99SBarry Smith } 146847c6ae99SBarry Smith 146947c6ae99SBarry Smith #undef __FUNCT__ 147047c6ae99SBarry Smith #define __FUNCT__ "DMSetFunction" 14717e833e3aSBarry Smith /*@C 147247c6ae99SBarry Smith DMSetFunction - sets a function to compute the right hand side vector entries for the KSP solver or nonlinear function for SNES 147347c6ae99SBarry Smith 147447c6ae99SBarry Smith Logically Collective on DM 147547c6ae99SBarry Smith 147647c6ae99SBarry Smith Input Parameter: 147747c6ae99SBarry Smith + dm - the DM object 147847c6ae99SBarry Smith - f - the function to compute (use PETSC_NULL to cancel a previous function that was set) 147947c6ae99SBarry Smith 148047c6ae99SBarry Smith Level: intermediate 148147c6ae99SBarry Smith 148247c6ae99SBarry Smith Notes: This sets both the function for function evaluations and the function used to compute Jacobians via finite differences if no Jacobian 148347c6ae99SBarry Smith computer is provided with DMSetJacobian(). Canceling cancels the function, but not the function used to compute the Jacobian. 148447c6ae99SBarry Smith 1485e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 148647c6ae99SBarry Smith DMSetJacobian() 148747c6ae99SBarry Smith 1488f07f9ceaSJed Brown @*/ 14897087cfbeSBarry Smith PetscErrorCode DMSetFunction(DM dm,PetscErrorCode (*f)(DM,Vec,Vec)) 149047c6ae99SBarry Smith { 149147c6ae99SBarry Smith PetscFunctionBegin; 1492171400e9SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 149347c6ae99SBarry Smith dm->ops->function = f; 149447c6ae99SBarry Smith if (f) { 149547c6ae99SBarry Smith dm->ops->functionj = f; 149647c6ae99SBarry Smith } 149747c6ae99SBarry Smith PetscFunctionReturn(0); 149847c6ae99SBarry Smith } 149947c6ae99SBarry Smith 150047c6ae99SBarry Smith #undef __FUNCT__ 150147c6ae99SBarry Smith #define __FUNCT__ "DMSetJacobian" 15027e833e3aSBarry Smith /*@C 150347c6ae99SBarry Smith DMSetJacobian - sets a function to compute the matrix entries for the KSP solver or Jacobian for SNES 150447c6ae99SBarry Smith 150547c6ae99SBarry Smith Logically Collective on DM 150647c6ae99SBarry Smith 150747c6ae99SBarry Smith Input Parameter: 150847c6ae99SBarry Smith + dm - the DM object to destroy 150947c6ae99SBarry Smith - f - the function to compute the matrix entries 151047c6ae99SBarry Smith 151147c6ae99SBarry Smith Level: intermediate 151247c6ae99SBarry Smith 1513e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 151447c6ae99SBarry Smith DMSetFunction() 151547c6ae99SBarry Smith 1516f07f9ceaSJed Brown @*/ 15177087cfbeSBarry Smith PetscErrorCode DMSetJacobian(DM dm,PetscErrorCode (*f)(DM,Vec,Mat,Mat,MatStructure*)) 151847c6ae99SBarry Smith { 151947c6ae99SBarry Smith PetscFunctionBegin; 1520171400e9SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 152147c6ae99SBarry Smith dm->ops->jacobian = f; 152247c6ae99SBarry Smith PetscFunctionReturn(0); 152347c6ae99SBarry Smith } 152447c6ae99SBarry Smith 152547c6ae99SBarry Smith #undef __FUNCT__ 152608da532bSDmitry Karpeev #define __FUNCT__ "DMSetVariableBounds" 152708da532bSDmitry Karpeev /*@C 152808da532bSDmitry Karpeev DMSetVariableBounds - sets a function to compute the the lower and upper bound vectors for SNESVI. 152908da532bSDmitry Karpeev 153008da532bSDmitry Karpeev Logically Collective on DM 153108da532bSDmitry Karpeev 153208da532bSDmitry Karpeev Input Parameter: 153308da532bSDmitry Karpeev + dm - the DM object 153408da532bSDmitry Karpeev - f - the function that computes variable bounds used by SNESVI (use PETSC_NULL to cancel a previous function that was set) 153508da532bSDmitry Karpeev 153608da532bSDmitry Karpeev Level: intermediate 153708da532bSDmitry Karpeev 1538e88d7f4bSDmitry Karpeev .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 153908da532bSDmitry Karpeev DMSetJacobian() 154008da532bSDmitry Karpeev 154108da532bSDmitry Karpeev @*/ 154208da532bSDmitry Karpeev PetscErrorCode DMSetVariableBounds(DM dm,PetscErrorCode (*f)(DM,Vec,Vec)) 154308da532bSDmitry Karpeev { 154408da532bSDmitry Karpeev PetscFunctionBegin; 154508da532bSDmitry Karpeev dm->ops->computevariablebounds = f; 154608da532bSDmitry Karpeev PetscFunctionReturn(0); 154708da532bSDmitry Karpeev } 154808da532bSDmitry Karpeev 154908da532bSDmitry Karpeev #undef __FUNCT__ 155008da532bSDmitry Karpeev #define __FUNCT__ "DMHasVariableBounds" 155108da532bSDmitry Karpeev /*@ 155208da532bSDmitry Karpeev DMHasVariableBounds - does the DM object have a variable bounds function? 155308da532bSDmitry Karpeev 155408da532bSDmitry Karpeev Not Collective 155508da532bSDmitry Karpeev 155608da532bSDmitry Karpeev Input Parameter: 155708da532bSDmitry Karpeev . dm - the DM object to destroy 155808da532bSDmitry Karpeev 155908da532bSDmitry Karpeev Output Parameter: 156008da532bSDmitry Karpeev . flg - PETSC_TRUE if the variable bounds function exists 156108da532bSDmitry Karpeev 156208da532bSDmitry Karpeev Level: developer 156308da532bSDmitry Karpeev 1564e88d7f4bSDmitry Karpeev .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian() 156508da532bSDmitry Karpeev 156608da532bSDmitry Karpeev @*/ 156708da532bSDmitry Karpeev PetscErrorCode DMHasVariableBounds(DM dm,PetscBool *flg) 156808da532bSDmitry Karpeev { 156908da532bSDmitry Karpeev PetscFunctionBegin; 157008da532bSDmitry Karpeev *flg = (dm->ops->computevariablebounds) ? PETSC_TRUE : PETSC_FALSE; 157108da532bSDmitry Karpeev PetscFunctionReturn(0); 157208da532bSDmitry Karpeev } 157308da532bSDmitry Karpeev 157408da532bSDmitry Karpeev #undef __FUNCT__ 157508da532bSDmitry Karpeev #define __FUNCT__ "DMComputeVariableBounds" 157608da532bSDmitry Karpeev /*@C 157708da532bSDmitry Karpeev DMComputeVariableBounds - compute variable bounds used by SNESVI. 157808da532bSDmitry Karpeev 157908da532bSDmitry Karpeev Logically Collective on DM 158008da532bSDmitry Karpeev 158108da532bSDmitry Karpeev Input Parameters: 158208da532bSDmitry Karpeev + dm - the DM object to destroy 158308da532bSDmitry Karpeev - x - current solution at which the bounds are computed 158408da532bSDmitry Karpeev 158508da532bSDmitry Karpeev Output parameters: 158608da532bSDmitry Karpeev + xl - lower bound 158708da532bSDmitry Karpeev - xu - upper bound 158808da532bSDmitry Karpeev 158908da532bSDmitry Karpeev Level: intermediate 159008da532bSDmitry Karpeev 1591e88d7f4bSDmitry Karpeev .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 159208da532bSDmitry Karpeev DMSetFunction(), DMSetVariableBounds() 159308da532bSDmitry Karpeev 159408da532bSDmitry Karpeev @*/ 159508da532bSDmitry Karpeev PetscErrorCode DMComputeVariableBounds(DM dm, Vec xl, Vec xu) 159608da532bSDmitry Karpeev { 159708da532bSDmitry Karpeev PetscErrorCode ierr; 159808da532bSDmitry Karpeev PetscFunctionBegin; 159908da532bSDmitry Karpeev PetscValidHeaderSpecific(xl,VEC_CLASSID,2); 160008da532bSDmitry Karpeev PetscValidHeaderSpecific(xu,VEC_CLASSID,2); 160108da532bSDmitry Karpeev if(dm->ops->computevariablebounds) { 160208da532bSDmitry Karpeev ierr = (*dm->ops->computevariablebounds)(dm, xl,xu); CHKERRQ(ierr); 160308da532bSDmitry Karpeev } 160408da532bSDmitry Karpeev else { 160508da532bSDmitry Karpeev ierr = VecSet(xl,SNES_VI_NINF); CHKERRQ(ierr); 160608da532bSDmitry Karpeev ierr = VecSet(xu,SNES_VI_INF); CHKERRQ(ierr); 160708da532bSDmitry Karpeev } 160808da532bSDmitry Karpeev PetscFunctionReturn(0); 160908da532bSDmitry Karpeev } 161008da532bSDmitry Karpeev 161108da532bSDmitry Karpeev #undef __FUNCT__ 161247c6ae99SBarry Smith #define __FUNCT__ "DMComputeInitialGuess" 161347c6ae99SBarry Smith /*@ 161447c6ae99SBarry Smith DMComputeInitialGuess - computes an initial guess vector entries for the KSP solvers 161547c6ae99SBarry Smith 161647c6ae99SBarry Smith Collective on DM 161747c6ae99SBarry Smith 161847c6ae99SBarry Smith Input Parameter: 161947c6ae99SBarry Smith + dm - the DM object to destroy 162047c6ae99SBarry Smith - x - the vector to hold the initial guess values 162147c6ae99SBarry Smith 162247c6ae99SBarry Smith Level: developer 162347c6ae99SBarry Smith 1624e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetRhs(), DMSetMat() 162547c6ae99SBarry Smith 162647c6ae99SBarry Smith @*/ 16277087cfbeSBarry Smith PetscErrorCode DMComputeInitialGuess(DM dm,Vec x) 162847c6ae99SBarry Smith { 162947c6ae99SBarry Smith PetscErrorCode ierr; 163047c6ae99SBarry Smith PetscFunctionBegin; 1631171400e9SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 163247c6ae99SBarry Smith if (!dm->ops->initialguess) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to provide function with DMSetInitialGuess()"); 163347c6ae99SBarry Smith ierr = (*dm->ops->initialguess)(dm,x);CHKERRQ(ierr); 163447c6ae99SBarry Smith PetscFunctionReturn(0); 163547c6ae99SBarry Smith } 163647c6ae99SBarry Smith 163747c6ae99SBarry Smith #undef __FUNCT__ 163847c6ae99SBarry Smith #define __FUNCT__ "DMHasInitialGuess" 163947c6ae99SBarry Smith /*@ 164047c6ae99SBarry Smith DMHasInitialGuess - does the DM object have an initial guess function 164147c6ae99SBarry Smith 164247c6ae99SBarry Smith Not Collective 164347c6ae99SBarry Smith 164447c6ae99SBarry Smith Input Parameter: 164547c6ae99SBarry Smith . dm - the DM object to destroy 164647c6ae99SBarry Smith 164747c6ae99SBarry Smith Output Parameter: 164847c6ae99SBarry Smith . flg - PETSC_TRUE if function exists 164947c6ae99SBarry Smith 165047c6ae99SBarry Smith Level: developer 165147c6ae99SBarry Smith 1652e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian() 165347c6ae99SBarry Smith 165447c6ae99SBarry Smith @*/ 16557087cfbeSBarry Smith PetscErrorCode DMHasInitialGuess(DM dm,PetscBool *flg) 165647c6ae99SBarry Smith { 165747c6ae99SBarry Smith PetscFunctionBegin; 1658171400e9SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 165947c6ae99SBarry Smith *flg = (dm->ops->initialguess) ? PETSC_TRUE : PETSC_FALSE; 166047c6ae99SBarry Smith PetscFunctionReturn(0); 166147c6ae99SBarry Smith } 166247c6ae99SBarry Smith 166347c6ae99SBarry Smith #undef __FUNCT__ 166447c6ae99SBarry Smith #define __FUNCT__ "DMHasFunction" 166547c6ae99SBarry Smith /*@ 166647c6ae99SBarry Smith DMHasFunction - does the DM object have a function 166747c6ae99SBarry Smith 166847c6ae99SBarry Smith Not Collective 166947c6ae99SBarry Smith 167047c6ae99SBarry Smith Input Parameter: 167147c6ae99SBarry Smith . dm - the DM object to destroy 167247c6ae99SBarry Smith 167347c6ae99SBarry Smith Output Parameter: 167447c6ae99SBarry Smith . flg - PETSC_TRUE if function exists 167547c6ae99SBarry Smith 167647c6ae99SBarry Smith Level: developer 167747c6ae99SBarry Smith 1678e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian() 167947c6ae99SBarry Smith 168047c6ae99SBarry Smith @*/ 16817087cfbeSBarry Smith PetscErrorCode DMHasFunction(DM dm,PetscBool *flg) 168247c6ae99SBarry Smith { 168347c6ae99SBarry Smith PetscFunctionBegin; 1684171400e9SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 168547c6ae99SBarry Smith *flg = (dm->ops->function) ? PETSC_TRUE : PETSC_FALSE; 168647c6ae99SBarry Smith PetscFunctionReturn(0); 168747c6ae99SBarry Smith } 168847c6ae99SBarry Smith 168947c6ae99SBarry Smith #undef __FUNCT__ 169047c6ae99SBarry Smith #define __FUNCT__ "DMHasJacobian" 169147c6ae99SBarry Smith /*@ 169247c6ae99SBarry Smith DMHasJacobian - does the DM object have a matrix function 169347c6ae99SBarry Smith 169447c6ae99SBarry Smith Not Collective 169547c6ae99SBarry Smith 169647c6ae99SBarry Smith Input Parameter: 169747c6ae99SBarry Smith . dm - the DM object to destroy 169847c6ae99SBarry Smith 169947c6ae99SBarry Smith Output Parameter: 170047c6ae99SBarry Smith . flg - PETSC_TRUE if function exists 170147c6ae99SBarry Smith 170247c6ae99SBarry Smith Level: developer 170347c6ae99SBarry Smith 1704e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian() 170547c6ae99SBarry Smith 170647c6ae99SBarry Smith @*/ 17077087cfbeSBarry Smith PetscErrorCode DMHasJacobian(DM dm,PetscBool *flg) 170847c6ae99SBarry Smith { 170947c6ae99SBarry Smith PetscFunctionBegin; 1710171400e9SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 171147c6ae99SBarry Smith *flg = (dm->ops->jacobian) ? PETSC_TRUE : PETSC_FALSE; 171247c6ae99SBarry Smith PetscFunctionReturn(0); 171347c6ae99SBarry Smith } 171447c6ae99SBarry Smith 171547c6ae99SBarry Smith #undef __FUNCT__ 171608da532bSDmitry Karpeev #define __FUNCT__ "DMSetVec" 1717748fac09SDmitry Karpeev /*@C 171808da532bSDmitry Karpeev DMSetVec - set the vector at which to compute residual, Jacobian and VI bounds, if the problem is nonlinear. 171908da532bSDmitry Karpeev 172008da532bSDmitry Karpeev Collective on DM 172108da532bSDmitry Karpeev 172208da532bSDmitry Karpeev Input Parameter: 172308da532bSDmitry Karpeev + dm - the DM object 1724e88d7f4bSDmitry Karpeev - x - location to compute residual and Jacobian, if PETSC_NULL is passed to those routines; will be PETSC_NULL for linear problems. 172508da532bSDmitry Karpeev 172608da532bSDmitry Karpeev Level: developer 172708da532bSDmitry Karpeev 1728e88d7f4bSDmitry Karpeev .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 172908da532bSDmitry Karpeev DMSetFunction(), DMSetJacobian(), DMSetVariableBounds() 173008da532bSDmitry Karpeev 173108da532bSDmitry Karpeev @*/ 173208da532bSDmitry Karpeev PetscErrorCode DMSetVec(DM dm,Vec x) 173308da532bSDmitry Karpeev { 173408da532bSDmitry Karpeev PetscErrorCode ierr; 173508da532bSDmitry Karpeev PetscFunctionBegin; 173608da532bSDmitry Karpeev if (x) { 173708da532bSDmitry Karpeev if (!dm->x) { 173808da532bSDmitry Karpeev ierr = DMCreateGlobalVector(dm,&dm->x);CHKERRQ(ierr); 173908da532bSDmitry Karpeev } 174008da532bSDmitry Karpeev ierr = VecCopy(x,dm->x);CHKERRQ(ierr); 174108da532bSDmitry Karpeev } 174208da532bSDmitry Karpeev else if(dm->x) { 174308da532bSDmitry Karpeev ierr = VecDestroy(&dm->x); CHKERRQ(ierr); 174408da532bSDmitry Karpeev } 174508da532bSDmitry Karpeev PetscFunctionReturn(0); 174608da532bSDmitry Karpeev } 174708da532bSDmitry Karpeev 174808da532bSDmitry Karpeev 174908da532bSDmitry Karpeev #undef __FUNCT__ 175047c6ae99SBarry Smith #define __FUNCT__ "DMComputeFunction" 175147c6ae99SBarry Smith /*@ 175247c6ae99SBarry Smith DMComputeFunction - computes the right hand side vector entries for the KSP solver or nonlinear function for SNES 175347c6ae99SBarry Smith 175447c6ae99SBarry Smith Collective on DM 175547c6ae99SBarry Smith 175647c6ae99SBarry Smith Input Parameter: 175747c6ae99SBarry Smith + dm - the DM object to destroy 175847c6ae99SBarry Smith . x - the location where the function is evaluationed, may be ignored for linear problems 175947c6ae99SBarry Smith - b - the vector to hold the right hand side entries 176047c6ae99SBarry Smith 176147c6ae99SBarry Smith Level: developer 176247c6ae99SBarry Smith 1763e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 176447c6ae99SBarry Smith DMSetJacobian() 176547c6ae99SBarry Smith 176647c6ae99SBarry Smith @*/ 17677087cfbeSBarry Smith PetscErrorCode DMComputeFunction(DM dm,Vec x,Vec b) 176847c6ae99SBarry Smith { 176947c6ae99SBarry Smith PetscErrorCode ierr; 177047c6ae99SBarry Smith PetscFunctionBegin; 1771171400e9SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 177247c6ae99SBarry Smith if (!dm->ops->function) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to provide function with DMSetFunction()"); 1773644e2e5bSBarry Smith PetscStackPush("DM user function"); 177447c6ae99SBarry Smith ierr = (*dm->ops->function)(dm,x,b);CHKERRQ(ierr); 1775644e2e5bSBarry Smith PetscStackPop; 177647c6ae99SBarry Smith PetscFunctionReturn(0); 177747c6ae99SBarry Smith } 177847c6ae99SBarry Smith 177947c6ae99SBarry Smith 178008da532bSDmitry Karpeev 178147c6ae99SBarry Smith #undef __FUNCT__ 178247c6ae99SBarry Smith #define __FUNCT__ "DMComputeJacobian" 178347c6ae99SBarry Smith /*@ 178447c6ae99SBarry Smith DMComputeJacobian - compute the matrix entries for the solver 178547c6ae99SBarry Smith 178647c6ae99SBarry Smith Collective on DM 178747c6ae99SBarry Smith 178847c6ae99SBarry Smith Input Parameter: 178947c6ae99SBarry Smith + dm - the DM object 1790cab2e9ccSBarry Smith . x - location to compute Jacobian at; will be PETSC_NULL for linear problems, for nonlinear problems if not provided then pulled from DM 179147c6ae99SBarry Smith . A - matrix that defines the operator for the linear solve 179247c6ae99SBarry Smith - B - the matrix used to construct the preconditioner 179347c6ae99SBarry Smith 179447c6ae99SBarry Smith Level: developer 179547c6ae99SBarry Smith 1796e727c939SJed Brown .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 179747c6ae99SBarry Smith DMSetFunction() 179847c6ae99SBarry Smith 179947c6ae99SBarry Smith @*/ 18007087cfbeSBarry Smith PetscErrorCode DMComputeJacobian(DM dm,Vec x,Mat A,Mat B,MatStructure *stflag) 180147c6ae99SBarry Smith { 180247c6ae99SBarry Smith PetscErrorCode ierr; 180347c6ae99SBarry Smith 180447c6ae99SBarry Smith PetscFunctionBegin; 1805171400e9SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 180647c6ae99SBarry Smith if (!dm->ops->jacobian) { 180747c6ae99SBarry Smith ISColoring coloring; 180847c6ae99SBarry Smith MatFDColoring fd; 180947c6ae99SBarry Smith 1810171400e9SBarry Smith ierr = DMCreateColoring(dm,dm->coloringtype,MATAIJ,&coloring);CHKERRQ(ierr); 181147c6ae99SBarry Smith ierr = MatFDColoringCreate(B,coloring,&fd);CHKERRQ(ierr); 1812fcfd50ebSBarry Smith ierr = ISColoringDestroy(&coloring);CHKERRQ(ierr); 181347c6ae99SBarry Smith ierr = MatFDColoringSetFunction(fd,(PetscErrorCode (*)(void))dm->ops->functionj,dm);CHKERRQ(ierr); 18140bdded8aSJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)fd,((PetscObject)dm)->prefix);CHKERRQ(ierr); 18150bdded8aSJed Brown ierr = MatFDColoringSetFromOptions(fd);CHKERRQ(ierr); 181671cd77b2SBarry Smith 181747c6ae99SBarry Smith dm->fd = fd; 181847c6ae99SBarry Smith dm->ops->jacobian = DMComputeJacobianDefault; 18192533e041SBarry Smith 182071cd77b2SBarry Smith /* don't know why this is needed */ 182171cd77b2SBarry Smith ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr); 182247c6ae99SBarry Smith } 182347c6ae99SBarry Smith if (!x) x = dm->x; 182447c6ae99SBarry Smith ierr = (*dm->ops->jacobian)(dm,x,A,B,stflag);CHKERRQ(ierr); 1825cab2e9ccSBarry Smith 182671cd77b2SBarry Smith /* if matrix depends on x; i.e. nonlinear problem, keep copy of input vector since needed by multigrid methods to generate coarse grid matrices */ 1827649052a6SBarry Smith if (x) { 1828cab2e9ccSBarry Smith if (!dm->x) { 182971cd77b2SBarry Smith ierr = DMCreateGlobalVector(dm,&dm->x);CHKERRQ(ierr); 1830cab2e9ccSBarry Smith } 1831cab2e9ccSBarry Smith ierr = VecCopy(x,dm->x);CHKERRQ(ierr); 1832649052a6SBarry Smith } 1833a8248277SBarry Smith if (A != B) { 1834a8248277SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1835a8248277SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1836a8248277SBarry Smith } 183747c6ae99SBarry Smith PetscFunctionReturn(0); 183847c6ae99SBarry Smith } 1839264ace61SBarry Smith 1840cab2e9ccSBarry Smith 1841264ace61SBarry Smith PetscFList DMList = PETSC_NULL; 1842264ace61SBarry Smith PetscBool DMRegisterAllCalled = PETSC_FALSE; 1843264ace61SBarry Smith 1844264ace61SBarry Smith #undef __FUNCT__ 1845264ace61SBarry Smith #define __FUNCT__ "DMSetType" 1846264ace61SBarry Smith /*@C 1847264ace61SBarry Smith DMSetType - Builds a DM, for a particular DM implementation. 1848264ace61SBarry Smith 1849264ace61SBarry Smith Collective on DM 1850264ace61SBarry Smith 1851264ace61SBarry Smith Input Parameters: 1852264ace61SBarry Smith + dm - The DM object 1853264ace61SBarry Smith - method - The name of the DM type 1854264ace61SBarry Smith 1855264ace61SBarry Smith Options Database Key: 1856264ace61SBarry Smith . -dm_type <type> - Sets the DM type; use -help for a list of available types 1857264ace61SBarry Smith 1858264ace61SBarry Smith Notes: 1859e1589f56SBarry Smith See "petsc/include/petscdm.h" for available DM types (for instance, DM1D, DM2D, or DM3D). 1860264ace61SBarry Smith 1861264ace61SBarry Smith Level: intermediate 1862264ace61SBarry Smith 1863264ace61SBarry Smith .keywords: DM, set, type 1864264ace61SBarry Smith .seealso: DMGetType(), DMCreate() 1865264ace61SBarry Smith @*/ 18667087cfbeSBarry Smith PetscErrorCode DMSetType(DM dm, const DMType method) 1867264ace61SBarry Smith { 1868264ace61SBarry Smith PetscErrorCode (*r)(DM); 1869264ace61SBarry Smith PetscBool match; 1870264ace61SBarry Smith PetscErrorCode ierr; 1871264ace61SBarry Smith 1872264ace61SBarry Smith PetscFunctionBegin; 1873264ace61SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID,1); 1874264ace61SBarry Smith ierr = PetscTypeCompare((PetscObject) dm, method, &match);CHKERRQ(ierr); 1875264ace61SBarry Smith if (match) PetscFunctionReturn(0); 1876264ace61SBarry Smith 1877264ace61SBarry Smith if (!DMRegisterAllCalled) {ierr = DMRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 18784b91b6eaSBarry Smith ierr = PetscFListFind(DMList, ((PetscObject)dm)->comm, method,PETSC_TRUE,(void (**)(void)) &r);CHKERRQ(ierr); 1879264ace61SBarry Smith if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown DM type: %s", method); 1880264ace61SBarry Smith 1881264ace61SBarry Smith if (dm->ops->destroy) { 1882264ace61SBarry Smith ierr = (*dm->ops->destroy)(dm);CHKERRQ(ierr); 1883b5c23020SJed Brown dm->ops->destroy = PETSC_NULL; 1884264ace61SBarry Smith } 1885264ace61SBarry Smith ierr = (*r)(dm);CHKERRQ(ierr); 1886264ace61SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)dm,method);CHKERRQ(ierr); 1887264ace61SBarry Smith PetscFunctionReturn(0); 1888264ace61SBarry Smith } 1889264ace61SBarry Smith 1890264ace61SBarry Smith #undef __FUNCT__ 1891264ace61SBarry Smith #define __FUNCT__ "DMGetType" 1892264ace61SBarry Smith /*@C 1893264ace61SBarry Smith DMGetType - Gets the DM type name (as a string) from the DM. 1894264ace61SBarry Smith 1895264ace61SBarry Smith Not Collective 1896264ace61SBarry Smith 1897264ace61SBarry Smith Input Parameter: 1898264ace61SBarry Smith . dm - The DM 1899264ace61SBarry Smith 1900264ace61SBarry Smith Output Parameter: 1901264ace61SBarry Smith . type - The DM type name 1902264ace61SBarry Smith 1903264ace61SBarry Smith Level: intermediate 1904264ace61SBarry Smith 1905264ace61SBarry Smith .keywords: DM, get, type, name 1906264ace61SBarry Smith .seealso: DMSetType(), DMCreate() 1907264ace61SBarry Smith @*/ 19087087cfbeSBarry Smith PetscErrorCode DMGetType(DM dm, const DMType *type) 1909264ace61SBarry Smith { 1910264ace61SBarry Smith PetscErrorCode ierr; 1911264ace61SBarry Smith 1912264ace61SBarry Smith PetscFunctionBegin; 1913264ace61SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID,1); 1914264ace61SBarry Smith PetscValidCharPointer(type,2); 1915264ace61SBarry Smith if (!DMRegisterAllCalled) { 1916264ace61SBarry Smith ierr = DMRegisterAll(PETSC_NULL);CHKERRQ(ierr); 1917264ace61SBarry Smith } 1918264ace61SBarry Smith *type = ((PetscObject)dm)->type_name; 1919264ace61SBarry Smith PetscFunctionReturn(0); 1920264ace61SBarry Smith } 1921264ace61SBarry Smith 192267a56275SMatthew G Knepley #undef __FUNCT__ 192367a56275SMatthew G Knepley #define __FUNCT__ "DMConvert" 192467a56275SMatthew G Knepley /*@C 192567a56275SMatthew G Knepley DMConvert - Converts a DM to another DM, either of the same or different type. 192667a56275SMatthew G Knepley 192767a56275SMatthew G Knepley Collective on DM 192867a56275SMatthew G Knepley 192967a56275SMatthew G Knepley Input Parameters: 193067a56275SMatthew G Knepley + dm - the DM 193167a56275SMatthew G Knepley - newtype - new DM type (use "same" for the same type) 193267a56275SMatthew G Knepley 193367a56275SMatthew G Knepley Output Parameter: 193467a56275SMatthew G Knepley . M - pointer to new DM 193567a56275SMatthew G Knepley 193667a56275SMatthew G Knepley Notes: 193767a56275SMatthew G Knepley Cannot be used to convert a sequential DM to parallel or parallel to sequential, 193867a56275SMatthew G Knepley the MPI communicator of the generated DM is always the same as the communicator 193967a56275SMatthew G Knepley of the input DM. 194067a56275SMatthew G Knepley 194167a56275SMatthew G Knepley Level: intermediate 194267a56275SMatthew G Knepley 194367a56275SMatthew G Knepley .seealso: DMCreate() 194467a56275SMatthew G Knepley @*/ 194567a56275SMatthew G Knepley PetscErrorCode DMConvert(DM dm, const DMType newtype, DM *M) 194667a56275SMatthew G Knepley { 194767a56275SMatthew G Knepley DM B; 194867a56275SMatthew G Knepley char convname[256]; 194967a56275SMatthew G Knepley PetscBool sametype, issame; 195067a56275SMatthew G Knepley PetscErrorCode ierr; 195167a56275SMatthew G Knepley 195267a56275SMatthew G Knepley PetscFunctionBegin; 195367a56275SMatthew G Knepley PetscValidHeaderSpecific(dm,DM_CLASSID,1); 195467a56275SMatthew G Knepley PetscValidType(dm,1); 195567a56275SMatthew G Knepley PetscValidPointer(M,3); 195667a56275SMatthew G Knepley ierr = PetscTypeCompare((PetscObject) dm, newtype, &sametype);CHKERRQ(ierr); 195767a56275SMatthew G Knepley ierr = PetscStrcmp(newtype, "same", &issame);CHKERRQ(ierr); 195867a56275SMatthew G Knepley { 195967a56275SMatthew G Knepley PetscErrorCode (*conv)(DM, const DMType, DM *) = PETSC_NULL; 196067a56275SMatthew G Knepley 196167a56275SMatthew G Knepley /* 196267a56275SMatthew G Knepley Order of precedence: 196367a56275SMatthew G Knepley 1) See if a specialized converter is known to the current DM. 196467a56275SMatthew G Knepley 2) See if a specialized converter is known to the desired DM class. 196567a56275SMatthew G Knepley 3) See if a good general converter is registered for the desired class 196667a56275SMatthew G Knepley 4) See if a good general converter is known for the current matrix. 196767a56275SMatthew G Knepley 5) Use a really basic converter. 196867a56275SMatthew G Knepley */ 196967a56275SMatthew G Knepley 197067a56275SMatthew G Knepley /* 1) See if a specialized converter is known to the current DM and the desired class */ 197167a56275SMatthew G Knepley ierr = PetscStrcpy(convname,"DMConvert_");CHKERRQ(ierr); 197267a56275SMatthew G Knepley ierr = PetscStrcat(convname,((PetscObject) dm)->type_name);CHKERRQ(ierr); 197367a56275SMatthew G Knepley ierr = PetscStrcat(convname,"_");CHKERRQ(ierr); 197467a56275SMatthew G Knepley ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr); 197567a56275SMatthew G Knepley ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr); 197667a56275SMatthew G Knepley ierr = PetscObjectQueryFunction((PetscObject)dm,convname,(void (**)(void))&conv);CHKERRQ(ierr); 197767a56275SMatthew G Knepley if (conv) goto foundconv; 197867a56275SMatthew G Knepley 197967a56275SMatthew G Knepley /* 2) See if a specialized converter is known to the desired DM class. */ 198067a56275SMatthew G Knepley ierr = DMCreate(((PetscObject) dm)->comm, &B);CHKERRQ(ierr); 198167a56275SMatthew G Knepley ierr = DMSetType(B, newtype);CHKERRQ(ierr); 198267a56275SMatthew G Knepley ierr = PetscStrcpy(convname,"DMConvert_");CHKERRQ(ierr); 198367a56275SMatthew G Knepley ierr = PetscStrcat(convname,((PetscObject) dm)->type_name);CHKERRQ(ierr); 198467a56275SMatthew G Knepley ierr = PetscStrcat(convname,"_");CHKERRQ(ierr); 198567a56275SMatthew G Knepley ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr); 198667a56275SMatthew G Knepley ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr); 198767a56275SMatthew G Knepley ierr = PetscObjectQueryFunction((PetscObject)B,convname,(void (**)(void))&conv);CHKERRQ(ierr); 198867a56275SMatthew G Knepley if (conv) { 1989fcfd50ebSBarry Smith ierr = DMDestroy(&B);CHKERRQ(ierr); 199067a56275SMatthew G Knepley goto foundconv; 199167a56275SMatthew G Knepley } 199267a56275SMatthew G Knepley 199367a56275SMatthew G Knepley #if 0 199467a56275SMatthew G Knepley /* 3) See if a good general converter is registered for the desired class */ 199567a56275SMatthew G Knepley conv = B->ops->convertfrom; 1996fcfd50ebSBarry Smith ierr = DMDestroy(&B);CHKERRQ(ierr); 199767a56275SMatthew G Knepley if (conv) goto foundconv; 199867a56275SMatthew G Knepley 199967a56275SMatthew G Knepley /* 4) See if a good general converter is known for the current matrix */ 200067a56275SMatthew G Knepley if (dm->ops->convert) { 200167a56275SMatthew G Knepley conv = dm->ops->convert; 200267a56275SMatthew G Knepley } 200367a56275SMatthew G Knepley if (conv) goto foundconv; 200467a56275SMatthew G Knepley #endif 200567a56275SMatthew G Knepley 200667a56275SMatthew G Knepley /* 5) Use a really basic converter. */ 200767a56275SMatthew G Knepley SETERRQ2(((PetscObject) dm)->comm, PETSC_ERR_SUP, "No conversion possible between DM types %s and %s", ((PetscObject) dm)->type_name, newtype); 200867a56275SMatthew G Knepley 200967a56275SMatthew G Knepley foundconv: 201067a56275SMatthew G Knepley ierr = PetscLogEventBegin(DM_Convert,dm,0,0,0);CHKERRQ(ierr); 201167a56275SMatthew G Knepley ierr = (*conv)(dm,newtype,M);CHKERRQ(ierr); 201267a56275SMatthew G Knepley ierr = PetscLogEventEnd(DM_Convert,dm,0,0,0);CHKERRQ(ierr); 201367a56275SMatthew G Knepley } 201467a56275SMatthew G Knepley ierr = PetscObjectStateIncrease((PetscObject) *M);CHKERRQ(ierr); 201567a56275SMatthew G Knepley PetscFunctionReturn(0); 201667a56275SMatthew G Knepley } 2017264ace61SBarry Smith 2018264ace61SBarry Smith /*--------------------------------------------------------------------------------------------------------------------*/ 2019264ace61SBarry Smith 2020264ace61SBarry Smith #undef __FUNCT__ 2021264ace61SBarry Smith #define __FUNCT__ "DMRegister" 2022264ace61SBarry Smith /*@C 2023264ace61SBarry Smith DMRegister - See DMRegisterDynamic() 2024264ace61SBarry Smith 2025264ace61SBarry Smith Level: advanced 2026264ace61SBarry Smith @*/ 20277087cfbeSBarry Smith PetscErrorCode DMRegister(const char sname[], const char path[], const char name[], PetscErrorCode (*function)(DM)) 2028264ace61SBarry Smith { 2029264ace61SBarry Smith char fullname[PETSC_MAX_PATH_LEN]; 2030264ace61SBarry Smith PetscErrorCode ierr; 2031264ace61SBarry Smith 2032264ace61SBarry Smith PetscFunctionBegin; 2033264ace61SBarry Smith ierr = PetscStrcpy(fullname, path);CHKERRQ(ierr); 2034264ace61SBarry Smith ierr = PetscStrcat(fullname, ":");CHKERRQ(ierr); 2035264ace61SBarry Smith ierr = PetscStrcat(fullname, name);CHKERRQ(ierr); 2036264ace61SBarry Smith ierr = PetscFListAdd(&DMList, sname, fullname, (void (*)(void)) function);CHKERRQ(ierr); 2037264ace61SBarry Smith PetscFunctionReturn(0); 2038264ace61SBarry Smith } 2039264ace61SBarry Smith 2040264ace61SBarry Smith 2041264ace61SBarry Smith /*--------------------------------------------------------------------------------------------------------------------*/ 2042264ace61SBarry Smith #undef __FUNCT__ 2043264ace61SBarry Smith #define __FUNCT__ "DMRegisterDestroy" 2044264ace61SBarry Smith /*@C 2045264ace61SBarry Smith DMRegisterDestroy - Frees the list of DM methods that were registered by DMRegister()/DMRegisterDynamic(). 2046264ace61SBarry Smith 2047264ace61SBarry Smith Not Collective 2048264ace61SBarry Smith 2049264ace61SBarry Smith Level: advanced 2050264ace61SBarry Smith 2051264ace61SBarry Smith .keywords: DM, register, destroy 2052264ace61SBarry Smith .seealso: DMRegister(), DMRegisterAll(), DMRegisterDynamic() 2053264ace61SBarry Smith @*/ 20547087cfbeSBarry Smith PetscErrorCode DMRegisterDestroy(void) 2055264ace61SBarry Smith { 2056264ace61SBarry Smith PetscErrorCode ierr; 2057264ace61SBarry Smith 2058264ace61SBarry Smith PetscFunctionBegin; 2059264ace61SBarry Smith ierr = PetscFListDestroy(&DMList);CHKERRQ(ierr); 2060264ace61SBarry Smith DMRegisterAllCalled = PETSC_FALSE; 2061264ace61SBarry Smith PetscFunctionReturn(0); 2062264ace61SBarry Smith } 206323f975d1SBarry Smith 206423f975d1SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 2065c6db04a5SJed Brown #include <mex.h> 206623f975d1SBarry Smith 20673014e516SBarry Smith typedef struct {char *funcname; char *jacname; mxArray *ctx;} DMMatlabContext; 206823f975d1SBarry Smith 206923f975d1SBarry Smith #undef __FUNCT__ 207023f975d1SBarry Smith #define __FUNCT__ "DMComputeFunction_Matlab" 207123f975d1SBarry Smith /* 207223f975d1SBarry Smith DMComputeFunction_Matlab - Calls the function that has been set with 207323f975d1SBarry Smith DMSetFunctionMatlab(). 207423f975d1SBarry Smith 207523f975d1SBarry Smith For linear problems x is null 207623f975d1SBarry Smith 207723f975d1SBarry Smith .seealso: DMSetFunction(), DMGetFunction() 207823f975d1SBarry Smith */ 20797087cfbeSBarry Smith PetscErrorCode DMComputeFunction_Matlab(DM dm,Vec x,Vec y) 208023f975d1SBarry Smith { 208123f975d1SBarry Smith PetscErrorCode ierr; 208223f975d1SBarry Smith DMMatlabContext *sctx; 208323f975d1SBarry Smith int nlhs = 1,nrhs = 4; 208423f975d1SBarry Smith mxArray *plhs[1],*prhs[4]; 208523f975d1SBarry Smith long long int lx = 0,ly = 0,ls = 0; 208623f975d1SBarry Smith 208723f975d1SBarry Smith PetscFunctionBegin; 208823f975d1SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 208923f975d1SBarry Smith PetscValidHeaderSpecific(y,VEC_CLASSID,3); 209023f975d1SBarry Smith PetscCheckSameComm(dm,1,y,3); 209123f975d1SBarry Smith 209223f975d1SBarry Smith /* call Matlab function in ctx with arguments x and y */ 20931b2093e4SBarry Smith ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr); 209423f975d1SBarry Smith ierr = PetscMemcpy(&ls,&dm,sizeof(dm));CHKERRQ(ierr); 209523f975d1SBarry Smith ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 20963014e516SBarry Smith ierr = PetscMemcpy(&ly,&y,sizeof(y));CHKERRQ(ierr); 209723f975d1SBarry Smith prhs[0] = mxCreateDoubleScalar((double)ls); 209823f975d1SBarry Smith prhs[1] = mxCreateDoubleScalar((double)lx); 209923f975d1SBarry Smith prhs[2] = mxCreateDoubleScalar((double)ly); 210023f975d1SBarry Smith prhs[3] = mxCreateString(sctx->funcname); 2101b807a863SBarry Smith ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscDMComputeFunctionInternal");CHKERRQ(ierr); 210223f975d1SBarry Smith ierr = mxGetScalar(plhs[0]);CHKERRQ(ierr); 210323f975d1SBarry Smith mxDestroyArray(prhs[0]); 210423f975d1SBarry Smith mxDestroyArray(prhs[1]); 210523f975d1SBarry Smith mxDestroyArray(prhs[2]); 210623f975d1SBarry Smith mxDestroyArray(prhs[3]); 210723f975d1SBarry Smith mxDestroyArray(plhs[0]); 210823f975d1SBarry Smith PetscFunctionReturn(0); 210923f975d1SBarry Smith } 211023f975d1SBarry Smith 211123f975d1SBarry Smith 211223f975d1SBarry Smith #undef __FUNCT__ 211323f975d1SBarry Smith #define __FUNCT__ "DMSetFunctionMatlab" 211423f975d1SBarry Smith /* 211523f975d1SBarry Smith DMSetFunctionMatlab - Sets the function evaluation routine 211623f975d1SBarry Smith 211723f975d1SBarry Smith */ 21187087cfbeSBarry Smith PetscErrorCode DMSetFunctionMatlab(DM dm,const char *func) 211923f975d1SBarry Smith { 212023f975d1SBarry Smith PetscErrorCode ierr; 212123f975d1SBarry Smith DMMatlabContext *sctx; 212223f975d1SBarry Smith 212323f975d1SBarry Smith PetscFunctionBegin; 2124171400e9SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 212523f975d1SBarry Smith /* currently sctx is memory bleed */ 21261b2093e4SBarry Smith ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr); 21273014e516SBarry Smith if (!sctx) { 212823f975d1SBarry Smith ierr = PetscMalloc(sizeof(DMMatlabContext),&sctx);CHKERRQ(ierr); 21293014e516SBarry Smith } 213023f975d1SBarry Smith ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr); 21311b2093e4SBarry Smith ierr = DMSetApplicationContext(dm,sctx);CHKERRQ(ierr); 213223f975d1SBarry Smith ierr = DMSetFunction(dm,DMComputeFunction_Matlab);CHKERRQ(ierr); 213323f975d1SBarry Smith PetscFunctionReturn(0); 213423f975d1SBarry Smith } 21353014e516SBarry Smith 21363014e516SBarry Smith #undef __FUNCT__ 21373014e516SBarry Smith #define __FUNCT__ "DMComputeJacobian_Matlab" 21383014e516SBarry Smith /* 21393014e516SBarry Smith DMComputeJacobian_Matlab - Calls the function that has been set with 21403014e516SBarry Smith DMSetJacobianMatlab(). 21413014e516SBarry Smith 21423014e516SBarry Smith For linear problems x is null 21433014e516SBarry Smith 21443014e516SBarry Smith .seealso: DMSetFunction(), DMGetFunction() 21453014e516SBarry Smith */ 21467087cfbeSBarry Smith PetscErrorCode DMComputeJacobian_Matlab(DM dm,Vec x,Mat A,Mat B,MatStructure *str) 21473014e516SBarry Smith { 21483014e516SBarry Smith PetscErrorCode ierr; 21493014e516SBarry Smith DMMatlabContext *sctx; 21503014e516SBarry Smith int nlhs = 2,nrhs = 5; 21513014e516SBarry Smith mxArray *plhs[2],*prhs[5]; 21523014e516SBarry Smith long long int lx = 0,lA = 0,lB = 0,ls = 0; 21533014e516SBarry Smith 21543014e516SBarry Smith PetscFunctionBegin; 21553014e516SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 21563014e516SBarry Smith PetscValidHeaderSpecific(A,MAT_CLASSID,3); 21573014e516SBarry Smith 2158e3c5b3baSBarry Smith /* call MATLAB function in ctx with arguments x, A, and B */ 21591b2093e4SBarry Smith ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr); 21603014e516SBarry Smith ierr = PetscMemcpy(&ls,&dm,sizeof(dm));CHKERRQ(ierr); 21613014e516SBarry Smith ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 21623014e516SBarry Smith ierr = PetscMemcpy(&lA,&A,sizeof(A));CHKERRQ(ierr); 21633014e516SBarry Smith ierr = PetscMemcpy(&lB,&B,sizeof(B));CHKERRQ(ierr); 21643014e516SBarry Smith prhs[0] = mxCreateDoubleScalar((double)ls); 21653014e516SBarry Smith prhs[1] = mxCreateDoubleScalar((double)lx); 21663014e516SBarry Smith prhs[2] = mxCreateDoubleScalar((double)lA); 21673014e516SBarry Smith prhs[3] = mxCreateDoubleScalar((double)lB); 21683014e516SBarry Smith prhs[4] = mxCreateString(sctx->jacname); 2169b807a863SBarry Smith ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscDMComputeJacobianInternal");CHKERRQ(ierr); 2170c980e822SBarry Smith *str = (MatStructure) mxGetScalar(plhs[0]); 2171c088a8dcSBarry Smith ierr = (PetscInt) mxGetScalar(plhs[1]);CHKERRQ(ierr); 21723014e516SBarry Smith mxDestroyArray(prhs[0]); 21733014e516SBarry Smith mxDestroyArray(prhs[1]); 21743014e516SBarry Smith mxDestroyArray(prhs[2]); 21753014e516SBarry Smith mxDestroyArray(prhs[3]); 21763014e516SBarry Smith mxDestroyArray(prhs[4]); 21773014e516SBarry Smith mxDestroyArray(plhs[0]); 21783014e516SBarry Smith mxDestroyArray(plhs[1]); 21793014e516SBarry Smith PetscFunctionReturn(0); 21803014e516SBarry Smith } 21813014e516SBarry Smith 21823014e516SBarry Smith 21833014e516SBarry Smith #undef __FUNCT__ 21843014e516SBarry Smith #define __FUNCT__ "DMSetJacobianMatlab" 21853014e516SBarry Smith /* 21863014e516SBarry Smith DMSetJacobianMatlab - Sets the Jacobian function evaluation routine 21873014e516SBarry Smith 21883014e516SBarry Smith */ 21897087cfbeSBarry Smith PetscErrorCode DMSetJacobianMatlab(DM dm,const char *func) 21903014e516SBarry Smith { 21913014e516SBarry Smith PetscErrorCode ierr; 21923014e516SBarry Smith DMMatlabContext *sctx; 21933014e516SBarry Smith 21943014e516SBarry Smith PetscFunctionBegin; 2195171400e9SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 21963014e516SBarry Smith /* currently sctx is memory bleed */ 21971b2093e4SBarry Smith ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr); 21983014e516SBarry Smith if (!sctx) { 21993014e516SBarry Smith ierr = PetscMalloc(sizeof(DMMatlabContext),&sctx);CHKERRQ(ierr); 22003014e516SBarry Smith } 22013014e516SBarry Smith ierr = PetscStrallocpy(func,&sctx->jacname);CHKERRQ(ierr); 22021b2093e4SBarry Smith ierr = DMSetApplicationContext(dm,sctx);CHKERRQ(ierr); 22033014e516SBarry Smith ierr = DMSetJacobian(dm,DMComputeJacobian_Matlab);CHKERRQ(ierr); 22043014e516SBarry Smith PetscFunctionReturn(0); 22053014e516SBarry Smith } 220623f975d1SBarry Smith #endif 2207b859378eSBarry Smith 2208b859378eSBarry Smith #undef __FUNCT__ 2209b859378eSBarry Smith #define __FUNCT__ "DMLoad" 2210b859378eSBarry Smith /*@C 2211b859378eSBarry Smith DMLoad - Loads a DM that has been stored in binary or HDF5 format 2212b859378eSBarry Smith with DMView(). 2213b859378eSBarry Smith 2214b859378eSBarry Smith Collective on PetscViewer 2215b859378eSBarry Smith 2216b859378eSBarry Smith Input Parameters: 2217b859378eSBarry Smith + newdm - the newly loaded DM, this needs to have been created with DMCreate() or 2218b859378eSBarry Smith some related function before a call to DMLoad(). 2219b859378eSBarry Smith - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or 2220b859378eSBarry Smith HDF5 file viewer, obtained from PetscViewerHDF5Open() 2221b859378eSBarry Smith 2222b859378eSBarry Smith Level: intermediate 2223b859378eSBarry Smith 2224b859378eSBarry Smith Notes: 2225b859378eSBarry Smith Defaults to the DM DA. 2226b859378eSBarry Smith 2227b859378eSBarry Smith Notes for advanced users: 2228b859378eSBarry Smith Most users should not need to know the details of the binary storage 2229b859378eSBarry Smith format, since DMLoad() and DMView() completely hide these details. 2230b859378eSBarry Smith But for anyone who's interested, the standard binary matrix storage 2231b859378eSBarry Smith format is 2232b859378eSBarry Smith .vb 2233b859378eSBarry Smith has not yet been determined 2234b859378eSBarry Smith .ve 2235b859378eSBarry Smith 2236b859378eSBarry Smith In addition, PETSc automatically does the byte swapping for 2237b859378eSBarry Smith machines that store the bytes reversed, e.g. DEC alpha, freebsd, 2238b859378eSBarry Smith linux, Windows and the paragon; thus if you write your own binary 2239b859378eSBarry Smith read/write routines you have to swap the bytes; see PetscBinaryRead() 2240b859378eSBarry Smith and PetscBinaryWrite() to see how this may be done. 2241b859378eSBarry Smith 2242b859378eSBarry Smith Concepts: vector^loading from file 2243b859378eSBarry Smith 2244b859378eSBarry Smith .seealso: PetscViewerBinaryOpen(), DMView(), MatLoad(), VecLoad() 2245b859378eSBarry Smith @*/ 2246b859378eSBarry Smith PetscErrorCode DMLoad(DM newdm, PetscViewer viewer) 2247b859378eSBarry Smith { 2248b859378eSBarry Smith PetscErrorCode ierr; 2249b859378eSBarry Smith 2250b859378eSBarry Smith PetscFunctionBegin; 2251b859378eSBarry Smith PetscValidHeaderSpecific(newdm,DM_CLASSID,1); 2252b859378eSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 2253b859378eSBarry Smith 2254b859378eSBarry Smith if (!((PetscObject)newdm)->type_name) { 2255b859378eSBarry Smith ierr = DMSetType(newdm, DMDA);CHKERRQ(ierr); 2256b859378eSBarry Smith } 2257b859378eSBarry Smith ierr = (*newdm->ops->load)(newdm,viewer);CHKERRQ(ierr); 2258b859378eSBarry Smith PetscFunctionReturn(0); 2259b859378eSBarry Smith } 2260b859378eSBarry Smith 22617da65231SMatthew G Knepley /******************************** FEM Support **********************************/ 22627da65231SMatthew G Knepley 22637da65231SMatthew G Knepley #undef __FUNCT__ 22647da65231SMatthew G Knepley #define __FUNCT__ "DMPrintCellVector" 22657da65231SMatthew G Knepley PetscErrorCode DMPrintCellVector(PetscInt c, const char name[], PetscInt len, const PetscScalar x[]) { 22661d47ebbbSSatish Balay PetscInt f; 22671b30c384SMatthew G Knepley PetscErrorCode ierr; 22681b30c384SMatthew G Knepley 22697da65231SMatthew G Knepley PetscFunctionBegin; 227074778d6cSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF, "Cell %D Element %s\n", c, name);CHKERRQ(ierr); 22711d47ebbbSSatish Balay for(f = 0; f < len; ++f) { 227274778d6cSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF, " | %G |\n", PetscRealPart(x[f]));CHKERRQ(ierr); 22737da65231SMatthew G Knepley } 22747da65231SMatthew G Knepley PetscFunctionReturn(0); 22757da65231SMatthew G Knepley } 22767da65231SMatthew G Knepley 22777da65231SMatthew G Knepley #undef __FUNCT__ 22787da65231SMatthew G Knepley #define __FUNCT__ "DMPrintCellMatrix" 22797da65231SMatthew G Knepley PetscErrorCode DMPrintCellMatrix(PetscInt c, const char name[], PetscInt rows, PetscInt cols, const PetscScalar A[]) { 22801b30c384SMatthew G Knepley PetscInt f, g; 22817da65231SMatthew G Knepley PetscErrorCode ierr; 22827da65231SMatthew G Knepley 22837da65231SMatthew G Knepley PetscFunctionBegin; 228474778d6cSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF, "Cell %D Element %s\n", c, name);CHKERRQ(ierr); 22851d47ebbbSSatish Balay for(f = 0; f < rows; ++f) { 228674778d6cSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF, " |");CHKERRQ(ierr); 22871d47ebbbSSatish Balay for(g = 0; g < cols; ++g) { 228874778d6cSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF, " % 9.5G", PetscRealPart(A[f*cols+g]));CHKERRQ(ierr); 22897da65231SMatthew G Knepley } 229074778d6cSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF, " |\n");CHKERRQ(ierr); 22917da65231SMatthew G Knepley } 22927da65231SMatthew G Knepley PetscFunctionReturn(0); 22937da65231SMatthew G Knepley } 2294e7c4fc90SDmitry Karpeev 2295