$extrastylesheet
00001 00002 #include "libmesh/petsc_macro.h" 00003 // This only works with petsc-3.3 and above. 00004 00005 #if !PETSC_VERSION_LESS_THAN(3,3,0) 00006 00007 // PETSc includes 00008 #include <petsc-private/dmimpl.h> 00009 00010 // Local Includes 00011 #include "libmesh/libmesh_common.h" 00012 #include "libmesh/nonlinear_implicit_system.h" 00013 #include "libmesh/petsc_nonlinear_solver.h" 00014 #include "libmesh/petsc_linear_solver.h" 00015 #include "libmesh/petsc_vector.h" 00016 #include "libmesh/petsc_matrix.h" 00017 #include "libmesh/petscdmlibmesh.h" 00018 #include "libmesh/dof_map.h" 00019 #include "libmesh/preconditioner.h" 00020 00021 00022 using namespace libMesh; 00023 00024 00025 #define DMLIBMESH_NO_DECOMPOSITION 0 00026 #define DMLIBMESH_FIELD_DECOMPOSITION 1 00027 #define DMLIBMESH_DOMAIN_DECOMPOSITION 2 00028 00029 #define DMLIBMESH_NO_EMBEDDING 0 00030 #define DMLIBMESH_FIELD_EMBEDDING 1 00031 #define DMLIBMESH_DOMAIN_EMBEDDING 2 00032 00033 struct DM_libMesh 00034 { 00035 NonlinearImplicitSystem* sys; 00036 std::map<std::string, unsigned int> *varids; 00037 std::map<unsigned int, std::string> *varnames; 00038 std::map<std::string, unsigned int> *blockids; 00039 std::map<unsigned int, std::string> *blocknames; 00040 unsigned int decomposition_type; 00041 std::vector<std::set<unsigned int> > *decomposition; 00042 unsigned int embedding_type; 00043 IS embedding; 00044 unsigned int vec_count; 00045 }; 00046 00047 struct DMVec_libMesh { 00048 std::string label; 00049 }; 00050 00051 #undef __FUNCT__ 00052 #define __FUNCT__ "DMlibMeshGetVec_Private" 00053 PetscErrorCode DMlibMeshGetVec_Private(DM /*dm*/, const char* /*name*/, Vec*) 00054 { 00055 PetscFunctionBegin; 00056 00057 PetscFunctionReturn(0); 00058 } 00059 00060 00061 00062 PetscErrorCode DMlibMeshSetUpName_Private(DM dm); 00063 00064 #undef __FUNCT__ 00065 #define __FUNCT__ "DMlibMeshSetSystem_libMesh" 00066 PetscErrorCode DMlibMeshSetSystem_libMesh(DM dm, NonlinearImplicitSystem& sys) 00067 { 00068 const Parallel::Communicator &comm(sys.comm()); 00069 00070 PetscErrorCode ierr; 00071 PetscFunctionBegin; 00072 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 00073 PetscBool islibmesh; 00074 ierr = PetscObjectTypeCompare((PetscObject)dm, DMLIBMESH,&islibmesh); 00075 if(!islibmesh) SETERRQ2(((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONG, "Got DM oftype %s, not of type %s", ((PetscObject)dm)->type_name, DMLIBMESH); 00076 00077 if(dm->setupcalled) SETERRQ(((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONGSTATE, "Cannot reset the libMesh system after DM has been set up."); 00078 DM_libMesh *dlm = (DM_libMesh *)(dm->data); 00079 dlm->sys =&sys; 00080 /* Initially populate the sets of active blockids and varids using all of the 00081 existing blocks/variables (only variables are supported at the moment). */ 00082 DofMap& dofmap = dlm->sys->get_dof_map(); 00083 dlm->varids->clear(); 00084 dlm->varnames->clear(); 00085 for(unsigned int v = 0; v < dofmap.n_variables(); ++v) { 00086 std::string vname = dofmap.variable(v).name(); 00087 dlm->varids->insert(std::pair<std::string,unsigned int>(vname,v)); 00088 dlm->varnames->insert(std::pair<unsigned int,std::string>(v,vname)); 00089 } 00090 const MeshBase& mesh = dlm->sys->get_mesh(); 00091 dlm->blockids->clear(); 00092 dlm->blocknames->clear(); 00093 std::set<subdomain_id_type> blocks; 00094 /* The following effectively is a verbatim copy of MeshBase::n_subdomains(). */ 00095 // This requires an inspection on every processor 00096 libmesh_parallel_only(mesh.comm()); 00097 MeshBase::const_element_iterator el = mesh.active_elements_begin(); 00098 const MeshBase::const_element_iterator end = mesh.active_elements_end(); 00099 for (; el!=end; ++el) 00100 blocks.insert((*el)->subdomain_id()); 00101 // Some subdomains may only live on other processors 00102 comm.set_union(blocks); 00103 00104 std::set<subdomain_id_type>::iterator bit = blocks.begin(); 00105 std::set<subdomain_id_type>::iterator bend = blocks.end(); 00106 if(bit == bend) SETERRQ(((PetscObject)dm)->comm, PETSC_ERR_PLIB, "No mesh blocks found."); 00107 00108 for(; bit != bend; ++bit) { 00109 subdomain_id_type bid = *bit; 00110 std::string bname = mesh.subdomain_name(bid); 00111 if(!bname.length()) { 00112 /* Block names are currently implemented for Exodus II meshes 00113 only, so we might have to make up our own block names and 00114 maintain our own mapping of block ids to names. 00115 */ 00116 std::ostringstream ss; 00117 ss << "dm" << bid; 00118 bname = ss.str(); 00119 } 00120 dlm->blockids->insert(std::pair<std::string,unsigned int>(bname,bid)); 00121 dlm->blocknames->insert(std::pair<unsigned int,std::string>(bid,bname)); 00122 } 00123 ierr = DMlibMeshSetUpName_Private(dm); CHKERRQ(ierr); 00124 PetscFunctionReturn(0); 00125 } 00126 00127 #undef __FUNCT__ 00128 #define __FUNCT__ "DMlibMeshGetSystem_libMesh" 00129 PetscErrorCode DMlibMeshGetSystem_libMesh(DM dm, NonlinearImplicitSystem*& sys) 00130 { 00131 PetscErrorCode ierr; 00132 00133 PetscFunctionBegin; 00134 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 00135 PetscBool islibmesh; 00136 ierr = PetscObjectTypeCompare((PetscObject)dm, DMLIBMESH,&islibmesh);CHKERRQ(ierr); 00137 if(!islibmesh) SETERRQ2(((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONG, "Got DM oftype %s, not of type %s", ((PetscObject)dm)->type_name, DMLIBMESH); 00138 DM_libMesh *dlm = (DM_libMesh *)(dm->data); 00139 sys = dlm->sys; 00140 PetscFunctionReturn(0); 00141 } 00142 00143 00144 #undef __FUNCT__ 00145 #define __FUNCT__ "DMlibMeshGetBlocks" 00146 PetscErrorCode DMlibMeshGetBlocks(DM dm, PetscInt *n, char*** blocknames) 00147 { 00148 PetscErrorCode ierr; 00149 PetscInt i; 00150 PetscFunctionBegin; 00151 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 00152 PetscBool islibmesh; 00153 ierr = PetscObjectTypeCompare((PetscObject)dm, DMLIBMESH,&islibmesh); 00154 if(!islibmesh) SETERRQ2(((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONG, "Got DM oftype %s, not of type %s", ((PetscObject)dm)->type_name, DMLIBMESH); 00155 DM_libMesh *dlm = (DM_libMesh *)(dm->data); 00156 PetscValidPointer(n,2); 00157 *n = dlm->blockids->size(); 00158 if(!blocknames) PetscFunctionReturn(0); 00159 ierr = PetscMalloc(*n*sizeof(char*), blocknames); CHKERRQ(ierr); 00160 i = 0; 00161 for(std::map<std::string, unsigned int>::const_iterator it = dlm->blockids->begin(); it != dlm->blockids->end(); ++it){ 00162 ierr = PetscStrallocpy(it->first.c_str(), *blocknames+i); CHKERRQ(ierr); 00163 ++i; 00164 } 00165 PetscFunctionReturn(0); 00166 } 00167 00168 #undef __FUNCT__ 00169 #define __FUNCT__ "DMlibMeshGetVariables" 00170 PetscErrorCode DMlibMeshGetVariables(DM dm, PetscInt *n, char*** varnames) 00171 { 00172 PetscErrorCode ierr; 00173 PetscFunctionBegin; 00174 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 00175 PetscBool islibmesh; 00176 PetscInt i; 00177 ierr = PetscObjectTypeCompare((PetscObject)dm, DMLIBMESH,&islibmesh); 00178 if(!islibmesh) SETERRQ2(((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONG, "Got DM oftype %s, not of type %s", ((PetscObject)dm)->type_name, DMLIBMESH); 00179 DM_libMesh *dlm = (DM_libMesh *)(dm->data); 00180 PetscValidPointer(n,2); 00181 *n = dlm->varids->size(); 00182 if(!varnames) PetscFunctionReturn(0); 00183 ierr = PetscMalloc(*n*sizeof(char*), varnames); CHKERRQ(ierr); 00184 i = 0; 00185 for(std::map<std::string, unsigned int>::const_iterator it = dlm->varids->begin(); it != dlm->varids->end(); ++it){ 00186 ierr = PetscStrallocpy(it->first.c_str(), *varnames+i); CHKERRQ(ierr); 00187 ++i; 00188 } 00189 PetscFunctionReturn(0); 00190 } 00191 00192 #undef __FUNCT__ 00193 #define __FUNCT__ "DMlibMeshSetUpName_Private" 00194 PetscErrorCode DMlibMeshSetUpName_Private(DM dm) 00195 { 00196 DM_libMesh* dlm = (DM_libMesh*)dm->data; 00197 PetscErrorCode ierr; 00198 PetscFunctionBegin; 00199 std::string name = dlm->sys->name(); 00200 std::map<unsigned int, std::string> *dnames = PETSC_NULL, *enames = PETSC_NULL; 00201 if(dlm->decomposition_type == DMLIBMESH_FIELD_DECOMPOSITION) { 00202 name += ":dec:var:"; 00203 dnames = dlm->varnames; 00204 } 00205 if(dlm->decomposition_type == DMLIBMESH_DOMAIN_DECOMPOSITION) { 00206 name += ":dec:block:"; 00207 dnames = dlm->blocknames; 00208 } 00209 if(dnames) { 00210 for(unsigned int d = 0; d < dlm->decomposition->size(); ++d) { 00211 for(std::set<unsigned int>::iterator dit = (*dlm->decomposition)[d].begin(); dit != (*dlm->decomposition)[d].end(); ++dit) { 00212 unsigned int id = *dit; 00213 if(dit != (*dlm->decomposition)[d].begin()) 00214 name += ","; 00215 name += (*dnames)[id]; 00216 } 00217 name += ";"; 00218 } 00219 } 00220 if(dlm->embedding_type == DMLIBMESH_FIELD_EMBEDDING) { 00221 name += ":emb:var:"; 00222 enames = dlm->varnames; 00223 } 00224 if(dlm->embedding_type == DMLIBMESH_DOMAIN_EMBEDDING) { 00225 name += ":emb:block:"; 00226 enames = dlm->blocknames; 00227 } 00228 if(enames) { 00229 for(std::map<unsigned int, std::string>::iterator eit = enames->begin(); eit != enames->end(); ++eit) { 00230 std::string ename = eit->second; 00231 if(eit != enames->begin()) 00232 name += ","; 00233 name += ename; 00234 } 00235 name += ";"; 00236 } 00237 ierr = PetscObjectSetName((PetscObject)dm, name.c_str()); CHKERRQ(ierr); 00238 PetscFunctionReturn(0); 00239 } 00240 00241 00242 #undef __FUNCT__ 00243 #define __FUNCT__ "DMCreateFieldDecomposition_libMesh" 00244 static PetscErrorCode DMCreateFieldDecomposition_libMesh(DM dm, PetscInt *len, char ***namelist, IS **islist, DM **dmlist) 00245 { 00246 PetscFunctionBegin; 00247 PetscErrorCode ierr; 00248 DM_libMesh *dlm = (DM_libMesh *)(dm->data); 00249 NonlinearImplicitSystem *sys = dlm->sys; 00250 IS emb; 00251 if(dlm->decomposition_type != DMLIBMESH_FIELD_DECOMPOSITION) PetscFunctionReturn(0); 00252 00253 *len = dlm->decomposition->size(); 00254 if(namelist) {ierr = PetscMalloc(*len*sizeof(char*), namelist); CHKERRQ(ierr);} 00255 if(islist) {ierr = PetscMalloc(*len*sizeof(IS), islist); CHKERRQ(ierr);} 00256 if(dmlist) {ierr = PetscMalloc(*len*sizeof(DM), dmlist); CHKERRQ(ierr);} 00257 DofMap& dofmap = dlm->sys->get_dof_map(); 00258 for(unsigned int d = 0; d < dlm->decomposition->size(); ++d) { 00259 std::set<numeric_index_type> dindices; 00260 std::string dname; 00261 std::map<std::string, unsigned int> dvarids; 00262 std::map<unsigned int, std::string> dvarnames; 00263 unsigned int dvcount = 0; 00264 for(std::set<unsigned int>::const_iterator dvit = (*dlm->decomposition)[d].begin(); dvit != (*dlm->decomposition)[d].end(); ++dvit){ 00265 unsigned int v = *dvit; 00266 std::string vname = (*dlm->varnames)[v]; 00267 dvarids.insert(std::pair<std::string, unsigned int>(vname,v)); 00268 dvarnames.insert(std::pair<unsigned int,std::string>(v,vname)); 00269 if(!dvcount) dname = vname; 00270 else dname += "_" + vname; 00271 ++dvcount; 00272 if(!islist) continue; 00273 /* Iterate only over this DM's blocks. */ 00274 for(std::map<std::string, unsigned int>::const_iterator bit = dlm->blockids->begin(); bit != dlm->blockids->end(); ++bit) { 00275 unsigned int b = bit->second; 00276 MeshBase::const_element_iterator el = sys->get_mesh().active_local_subdomain_elements_begin(b); 00277 MeshBase::const_element_iterator end_el = sys->get_mesh().active_local_subdomain_elements_end(b); 00278 for ( ; el != end_el; ++el) { 00279 const Elem* elem = *el; 00280 //unsigned int e_subdomain = elem->subdomain_id(); 00281 std::vector<numeric_index_type> evindices; 00282 // Get the degree of freedom indices for the given variable off the current element. 00283 dofmap.dof_indices(elem, evindices, v); 00284 for(unsigned int i = 0; i < evindices.size(); ++i) { 00285 numeric_index_type dof = evindices[i]; 00286 if(dof >= dofmap.first_dof() && dof < dofmap.end_dof()) /* might want to use variable_first/last_local_dof instead */ 00287 dindices.insert(dof); 00288 } 00289 } 00290 } 00291 } 00292 if(namelist) { 00293 ierr = PetscStrallocpy(dname.c_str(),(*namelist)+d); CHKERRQ(ierr); 00294 } 00295 if(islist) { 00296 IS dis; 00297 PetscInt *darray; 00298 ierr = PetscMalloc(sizeof(PetscInt)*dindices.size(), &darray); CHKERRQ(ierr); 00299 numeric_index_type i = 0; 00300 for(std::set<numeric_index_type>::const_iterator it = dindices.begin(); it != dindices.end(); ++it) { 00301 darray[i] = *it; 00302 ++i; 00303 } 00304 ierr = ISCreateGeneral(((PetscObject)dm)->comm, dindices.size(),darray, PETSC_OWN_POINTER, &dis); CHKERRQ(ierr); 00305 if(dlm->embedding) { 00306 /* Create a relative embedding into the parent's index space. */ 00307 #if PETSC_RELEASE_LESS_THAN(3,3,1) 00308 ierr = ISMapFactorRight(dis,dlm->embedding, PETSC_TRUE, &emb); CHKERRQ(ierr); 00309 #else 00310 ierr = ISEmbed(dis,dlm->embedding, PETSC_TRUE, &emb); CHKERRQ(ierr); 00311 #endif 00312 PetscInt elen, dlen; 00313 ierr = ISGetLocalSize(emb, &elen); CHKERRQ(ierr); 00314 ierr = ISGetLocalSize(dis, &dlen); CHKERRQ(ierr); 00315 if(elen != dlen) SETERRQ1(((PetscObject)dm)->comm, PETSC_ERR_PLIB, "Failed to embed subdomain %D", d); 00316 ierr = ISDestroy(&dis); CHKERRQ(ierr); 00317 dis = emb; 00318 } 00319 else { 00320 emb = dis; 00321 } 00322 (*islist)[d] = dis; 00323 } 00324 if(dmlist) { 00325 DM ddm; 00326 ierr = DMCreate(((PetscObject)dm)->comm, &ddm); CHKERRQ(ierr); 00327 ierr = DMSetType(ddm, DMLIBMESH); CHKERRQ(ierr); 00328 DM_libMesh *ddlm = (DM_libMesh*)(ddm->data); 00329 ddlm->sys = dlm->sys; 00330 /* copy over the block ids and names */ 00331 *ddlm->blockids = *dlm->blockids; 00332 *ddlm->blocknames = *dlm->blocknames; 00333 /* set the vars from the d-th part of the decomposition. */ 00334 *ddlm->varids = dvarids; 00335 *ddlm->varnames = dvarnames; 00336 ierr = PetscObjectReference((PetscObject)emb); CHKERRQ(ierr); 00337 ddlm->embedding = emb; 00338 ddlm->embedding_type = DMLIBMESH_FIELD_EMBEDDING; 00339 00340 ierr = DMlibMeshSetUpName_Private(ddm); CHKERRQ(ierr); 00341 ierr = DMSetFromOptions(ddm); CHKERRQ(ierr); 00342 (*dmlist)[d] = ddm; 00343 } 00344 } 00345 PetscFunctionReturn(0); 00346 } 00347 00348 #undef __FUNCT__ 00349 #define __FUNCT__ "DMCreateDomainDecomposition_libMesh" 00350 static PetscErrorCode DMCreateDomainDecomposition_libMesh(DM dm, PetscInt *len, char ***namelist, IS **innerislist, IS **outerislist, DM **dmlist) 00351 { 00352 PetscFunctionBegin; 00353 PetscErrorCode ierr; 00354 DM_libMesh *dlm = (DM_libMesh *)(dm->data); 00355 NonlinearImplicitSystem *sys = dlm->sys; 00356 IS emb; 00357 if(dlm->decomposition_type != DMLIBMESH_DOMAIN_DECOMPOSITION) PetscFunctionReturn(0); 00358 *len = dlm->decomposition->size(); 00359 if(namelist) {ierr = PetscMalloc(*len*sizeof(char*), namelist); CHKERRQ(ierr);} 00360 if(innerislist) {ierr = PetscMalloc(*len*sizeof(IS), innerislist); CHKERRQ(ierr);} 00361 if(outerislist) *outerislist = PETSC_NULL; /* FIX: allow mesh-based overlap. */ 00362 if(dmlist) {ierr = PetscMalloc(*len*sizeof(DM), dmlist); CHKERRQ(ierr);} 00363 for(unsigned int d = 0; d < dlm->decomposition->size(); ++d) { 00364 std::set<numeric_index_type> dindices; 00365 std::string dname; 00366 std::map<std::string, unsigned int> dblockids; 00367 std::map<unsigned int,std::string> dblocknames; 00368 unsigned int dbcount = 0; 00369 for(std::set<unsigned int>::const_iterator bit = (*dlm->decomposition)[d].begin(); bit != (*dlm->decomposition)[d].end(); ++bit){ 00370 unsigned int b = *bit; 00371 std::string bname = (*dlm->blocknames)[b]; 00372 dblockids.insert(std::pair<std::string, unsigned int>(bname,b)); 00373 dblocknames.insert(std::pair<unsigned int,std::string>(b,bname)); 00374 if(!dbcount) dname = bname; 00375 else dname += "_" + bname; 00376 ++dbcount; 00377 if(!innerislist) continue; 00378 MeshBase::const_element_iterator el = sys->get_mesh().active_local_subdomain_elements_begin(b); 00379 const MeshBase::const_element_iterator end_el = sys->get_mesh().active_local_subdomain_elements_end(b); 00380 for ( ; el != end_el; ++el) { 00381 const Elem* elem = *el; 00382 std::vector<numeric_index_type> evindices; 00383 /* Iterate only over this DM's variables. */ 00384 for(std::map<std::string, unsigned int>::const_iterator vit = dlm->varids->begin(); vit != dlm->varids->end(); ++vit) { 00385 unsigned int v = vit->second; 00386 // Get the degree of freedom indices for the given variable off the current element. 00387 sys->get_dof_map().dof_indices(elem, evindices, v); 00388 for(unsigned int i = 0; i < evindices.size(); ++i) { 00389 numeric_index_type dof = evindices[i]; 00390 if(dof >= sys->get_dof_map().first_dof() && dof < sys->get_dof_map().end_dof()) /* might want to use variable_first/last_local_dof instead */ 00391 dindices.insert(dof); 00392 } 00393 } 00394 } 00395 } 00396 if(namelist) { 00397 ierr = PetscStrallocpy(dname.c_str(),(*namelist)+d); CHKERRQ(ierr); 00398 } 00399 if(innerislist) { 00400 PetscInt *darray; 00401 IS dis; 00402 ierr = PetscMalloc(sizeof(PetscInt)*dindices.size(), &darray); CHKERRQ(ierr); 00403 numeric_index_type i = 0; 00404 for(std::set<numeric_index_type>::const_iterator it = dindices.begin(); it != dindices.end(); ++it) { 00405 darray[i] = *it; 00406 ++i; 00407 } 00408 ierr = ISCreateGeneral(((PetscObject)dm)->comm, dindices.size(),darray, PETSC_OWN_POINTER, &dis); CHKERRQ(ierr); 00409 if(dlm->embedding) { 00410 /* Create a relative embedding into the parent's index space. */ 00411 #if PETSC_RELEASE_LESS_THAN(3,3,1) 00412 ierr = ISMapFactorRight(dis,dlm->embedding, PETSC_TRUE, &emb); CHKERRQ(ierr); 00413 #else 00414 ierr = ISEmbed(dis,dlm->embedding, PETSC_TRUE, &emb); CHKERRQ(ierr); 00415 #endif 00416 PetscInt elen, dlen; 00417 ierr = ISGetLocalSize(emb, &elen); CHKERRQ(ierr); 00418 ierr = ISGetLocalSize(dis, &dlen); CHKERRQ(ierr); 00419 if(elen != dlen) SETERRQ1(((PetscObject)dm)->comm, PETSC_ERR_PLIB, "Failed to embed field %D", d); 00420 ierr = ISDestroy(&dis); CHKERRQ(ierr); 00421 dis = emb; 00422 } 00423 else { 00424 emb = dis; 00425 } 00426 if(innerislist) { 00427 ierr = PetscObjectReference((PetscObject)dis); CHKERRQ(ierr); 00428 (*innerislist)[d] = dis; 00429 } 00430 ierr = ISDestroy(&dis); CHKERRQ(ierr); 00431 } 00432 if(dmlist) { 00433 DM ddm; 00434 ierr = DMCreate(((PetscObject)dm)->comm, &ddm); CHKERRQ(ierr); 00435 ierr = DMSetType(ddm, DMLIBMESH); CHKERRQ(ierr); 00436 DM_libMesh *ddlm = (DM_libMesh*)(ddm->data); 00437 ddlm->sys = dlm->sys; 00438 /* copy over the varids and varnames */ 00439 *ddlm->varids = *dlm->varids; 00440 *ddlm->varnames = *dlm->varnames; 00441 /* set the blocks from the d-th part of the decomposition. */ 00442 *ddlm->blockids = dblockids; 00443 *ddlm->blocknames = dblocknames; 00444 ierr = PetscObjectReference((PetscObject)emb); CHKERRQ(ierr); 00445 ddlm->embedding = emb; 00446 ddlm->embedding_type = DMLIBMESH_DOMAIN_EMBEDDING; 00447 00448 ierr = DMlibMeshSetUpName_Private(ddm); CHKERRQ(ierr); 00449 ierr = DMSetFromOptions(ddm); CHKERRQ(ierr); 00450 (*dmlist)[d] = ddm; 00451 } 00452 } 00453 PetscFunctionReturn(0); 00454 } 00455 00456 00457 #undef __FUNCT__ 00458 #define __FUNCT__ "DMlibMeshCreateFieldDecompositionDM" 00459 PetscErrorCode DMlibMeshCreateFieldDecompositionDM(DM dm, PetscInt dnumber, PetscInt* dsizes, char*** dvarlists, DM* ddm) 00460 { 00461 PetscErrorCode ierr; 00462 PetscBool islibmesh; 00463 PetscFunctionBegin; 00464 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 00465 ierr = PetscObjectTypeCompare((PetscObject)dm, DMLIBMESH,&islibmesh); 00466 if(!islibmesh) SETERRQ2(((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONG, "Got DM oftype %s, not of type %s", ((PetscObject)dm)->type_name, DMLIBMESH); 00467 if(dnumber < 0) SETERRQ1(((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONG, "Negative number %D of decomposition parts", dnumber); 00468 PetscValidPointer(ddm,5); 00469 DM_libMesh *dlm = (DM_libMesh *)(dm->data); 00470 ierr = DMCreate(((PetscObject)dm)->comm, ddm); CHKERRQ(ierr); 00471 ierr = DMSetType(*ddm, DMLIBMESH); CHKERRQ(ierr); 00472 DM_libMesh *ddlm = (DM_libMesh *)((*ddm)->data); 00473 ddlm->sys = dlm->sys; 00474 ddlm->varids = dlm->varids; 00475 ddlm->varnames = dlm->varnames; 00476 ddlm->blockids = dlm->blockids; 00477 ddlm->blocknames = dlm->blocknames; 00478 ddlm->decomposition = new(std::vector<std::set<unsigned int> >); 00479 ddlm->decomposition_type = DMLIBMESH_FIELD_DECOMPOSITION; 00480 if(dnumber) { 00481 for(PetscInt d = 0; d < dnumber; ++d) { 00482 if(dsizes[d] < 0) SETERRQ2(((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONG, "Negative size %D of decomposition part %D", dsizes[d],d); 00483 ddlm->decomposition->push_back(std::set<unsigned int>()); 00484 for(PetscInt v = 0; v < dsizes[d]; ++v) { 00485 std::string vname(dvarlists[d][v]); 00486 std::map<std::string, unsigned int>::const_iterator vit = dlm->varids->find(vname); 00487 if(vit == dlm->varids->end()) 00488 SETERRQ3(((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONG, "Variable %D on the %D-th list with name %s is not owned by this DM", v, d, dvarlists[d][v]); 00489 unsigned int vid = vit->second; 00490 (*ddlm->decomposition)[d].insert(vid); 00491 } 00492 } 00493 } 00494 else { /* Empty splits indicate default: split all variables with one per split. */ 00495 PetscInt d = 0; 00496 for(std::map<std::string, unsigned int>::const_iterator vit = ddlm->varids->begin(); vit != ddlm->varids->end(); ++vit) { 00497 ddlm->decomposition->push_back(std::set<unsigned int>()); 00498 unsigned int vid = vit->second; 00499 std::string vname = vit->first; 00500 (*ddlm->decomposition)[d].insert(vid); 00501 ++d; 00502 } 00503 } 00504 ierr = DMlibMeshSetUpName_Private(*ddm); CHKERRQ(ierr); 00505 ierr = DMSetFromOptions(*ddm); CHKERRQ(ierr); 00506 ierr = DMSetUp(*ddm); CHKERRQ(ierr); 00507 PetscFunctionReturn(0); 00508 } 00509 00510 #undef __FUNCT__ 00511 #define __FUNCT__ "DMlibMeshCreateDomainDecompositionDM" 00512 PetscErrorCode DMlibMeshCreateDomainDecompositionDM(DM dm, PetscInt dnumber, PetscInt* dsizes, char*** dblocklists, DM* ddm) 00513 { 00514 PetscErrorCode ierr; 00515 PetscBool islibmesh; 00516 PetscFunctionBegin; 00517 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 00518 ierr = PetscObjectTypeCompare((PetscObject)dm, DMLIBMESH,&islibmesh); 00519 if(!islibmesh) SETERRQ2(((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONG, "Got DM oftype %s, not of type %s", ((PetscObject)dm)->type_name, DMLIBMESH); 00520 if(dnumber < 0) SETERRQ1(((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONG, "Negative number %D of decomposition parts", dnumber); 00521 PetscValidPointer(ddm,5); 00522 DM_libMesh *dlm = (DM_libMesh *)(dm->data); 00523 ierr = DMCreate(((PetscObject)dm)->comm, ddm); CHKERRQ(ierr); 00524 ierr = DMSetType(*ddm, DMLIBMESH); CHKERRQ(ierr); 00525 DM_libMesh *ddlm = (DM_libMesh *)((*ddm)->data); 00526 ddlm->sys = dlm->sys; 00527 ddlm->varids = dlm->varids; 00528 ddlm->varnames = dlm->varnames; 00529 ddlm->blockids = dlm->blockids; 00530 ddlm->blocknames = dlm->blocknames; 00531 ddlm->decomposition = new(std::vector<std::set<unsigned int> >); 00532 ddlm->decomposition_type = DMLIBMESH_DOMAIN_DECOMPOSITION; 00533 if(dnumber) { 00534 for(PetscInt d = 0; d < dnumber; ++d) { 00535 if(dsizes[d] < 0) SETERRQ2(((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONG, "Negative size %D of decomposition part %D", dsizes[d],d); 00536 ddlm->decomposition->push_back(std::set<unsigned int>()); 00537 for(PetscInt b = 0; b < dsizes[d]; ++b) { 00538 std::string bname(dblocklists[d][b]); 00539 std::map<std::string, unsigned int>::const_iterator bit = dlm->blockids->find(bname); 00540 if(bit == dlm->blockids->end()) 00541 SETERRQ3(((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONG, "Block %D on the %D-th list with name %s is not owned by this DM", b, d, dblocklists[d][b]); 00542 unsigned int bid = bit->second; 00543 (*ddlm->decomposition)[d].insert(bid); 00544 } 00545 } 00546 } 00547 else { /* Empty splits indicate default: split all blocks with one per split. */ 00548 PetscInt d = 0; 00549 for(std::map<std::string, unsigned int>::const_iterator bit = ddlm->blockids->begin(); bit != ddlm->blockids->end(); ++bit) { 00550 ddlm->decomposition->push_back(std::set<unsigned int>()); 00551 unsigned int bid = bit->second; 00552 std::string bname = bit->first; 00553 (*ddlm->decomposition)[d].insert(bid); 00554 ++d; 00555 } 00556 } 00557 ierr = DMlibMeshSetUpName_Private(*ddm); CHKERRQ(ierr); 00558 ierr = DMSetFromOptions(*ddm); CHKERRQ(ierr); 00559 ierr = DMSetUp(*ddm); CHKERRQ(ierr); 00560 PetscFunctionReturn(0); 00561 } 00562 00563 struct token { 00564 const char* s; 00565 struct token *next; 00566 }; 00567 00568 00569 // The following functions are only needed for older PETScs. 00570 #if PETSC_RELEASE_LESS_THAN(3,3,1) 00571 00572 #undef __FUNCT__ 00573 #define __FUNCT__ "DMlibMeshParseDecompositionDescriptor_Private" 00574 static PetscErrorCode DMlibMeshParseDecompositionDescriptor_Private(DM dm, const char* ddesc, PetscInt* dtype, PetscInt* dcount, PetscInt** dsizes, char ****dlists) 00575 { 00576 PetscFunctionBegin; 00577 PetscErrorCode ierr; 00578 PetscBool eq; 00579 char *s0, *s, *ss; 00580 struct token *llfirst = PETSC_NULL, *lllast = PETSC_NULL, *tok; 00581 PetscInt stcount = 0, brcount = 0, d, i; 00582 size_t len0, count; 00583 00584 /* 00585 Parse the decomposition descriptor. 00586 Decomposition names could be of one of two forms: 00587 var:v1,v2;v3,v4;v4,v5; 00588 block:b1,b2;b3,b4;b4,b5; 00589 resulting in an overlapping decomposition that groups 00590 variables (v1,v2), (v3,v4), (v4,v5) or 00591 blocks (b1,b2), (b3,b4), (b4,b5). 00592 */ 00593 /* Copy the descriptor so that we can manipulate it in place. */ 00594 ierr = PetscStrallocpy(ddesc,&s0); CHKERRQ(ierr); 00595 ierr = PetscStrlen(s0, &len0) ; CHKERRQ(ierr); 00596 ierr = PetscStrstr(s0,":",&ss); CHKERRQ(ierr); 00597 if(!ss) { 00598 ss = s0+len0; 00599 } 00600 else { 00601 *ss = 0; 00602 } 00603 ierr = PetscStrcmp(s0,"var",&eq); CHKERRQ(ierr); 00604 if(eq) { 00605 *dtype=DMLIBMESH_FIELD_DECOMPOSITION; 00606 } 00607 else { 00608 ierr = PetscStrcmp(s0,"block",&eq);CHKERRQ(ierr); 00609 if(!eq) 00610 SETERRQ1(((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONG, "Could not determine decomposition type from descriptor: %s\n", ddesc); CHKERRQ(ierr); 00611 *dtype=DMLIBMESH_DOMAIN_DECOMPOSITION; 00612 } 00613 ierr = PetscStrlen(s0,&count); CHKERRQ(ierr); 00614 while(count < len0) { 00615 struct token *st, *br; 00616 ++ss; ++count; 00617 s=ss; 00618 while(*ss && *ss != ',' && *ss != ';') { 00619 ++ss; ++count; 00620 } 00621 st = PETSC_NULL; br = PETSC_NULL; 00622 if(*ss) { 00623 /* 00624 Found a separator, or a break. 00625 Add an appropriate token to the list. 00626 A token separator ',' produces no token. 00627 */ 00628 if(*ss == ';') { 00629 /* Create a break token: a token with a null string. */ 00630 #if PETSC_RELEASE_LESS_THAN(3,5,0) 00631 ierr = PetscNew(struct token,&br);CHKERRQ(ierr); 00632 #else 00633 ierr = PetscNew(&br);CHKERRQ(ierr); 00634 #endif 00635 } 00636 *ss = 0; 00637 if(s != ss) { 00638 /* A nonempty string. */ 00639 #if PETSC_RELEASE_LESS_THAN(3,5,0) 00640 ierr = PetscNew(struct token, &st);CHKERRQ(ierr); 00641 #else 00642 ierr = PetscNew(&st);CHKERRQ(ierr); 00643 #endif 00644 st->s = s; /* The string will be properly copied below. */ 00645 } 00646 /* Add the new tokens to the list. */ 00647 if(st) { 00648 if(!lllast) { 00649 llfirst = lllast = st; 00650 } 00651 else { 00652 lllast->next = st; lllast = st; 00653 } 00654 } 00655 if(br) { 00656 if(!lllast) { 00657 llfirst = lllast = br; 00658 } 00659 else { 00660 lllast->next = br; lllast = br; 00661 } 00662 } 00663 } 00664 } 00665 /* The result of parsing is in the linked list ll. */ 00666 /* Count up the strings and the breaks. */ 00667 tok = llfirst; 00668 while(tok) { 00669 if(tok->s) 00670 ++stcount; 00671 else 00672 ++brcount; 00673 tok = tok->next; 00674 } 00675 /* Allocate the space for the output. */ 00676 *dcount = brcount; 00677 ierr = PetscMalloc(*dcount*sizeof(PetscInt), dsizes); CHKERRQ(ierr); 00678 ierr = PetscMalloc(*dcount*sizeof(char**), dlists); CHKERRQ(ierr); 00679 for(d = 0; d < *dcount; ++d) (*dsizes)[d] = 0; 00680 tok = llfirst; d = 0; 00681 while(tok) { 00682 if(tok->s) 00683 ++(*dsizes)[d]; 00684 else 00685 ++d; 00686 tok = tok->next; 00687 } 00688 for(d = 0; d < *dcount; ++d) { 00689 ierr = PetscMalloc(sizeof(char**)*(*dsizes)[d], (*dlists)+d); CHKERRQ(ierr); 00690 } 00691 /* Now copy strings and destroy tokens. */ 00692 tok = llfirst; d = 0; i = 0; 00693 while(tok) { 00694 if(tok->s) { 00695 ierr = PetscStrallocpy(tok->s, (*dlists)[d]+i); CHKERRQ(ierr); 00696 ++i; 00697 } 00698 else { 00699 ++d; 00700 i = 0; 00701 } 00702 llfirst = tok; 00703 tok = tok->next; 00704 ierr = PetscFree(llfirst); CHKERRQ(ierr); 00705 } 00706 /* Deallocate workspace. */ 00707 ierr = PetscFree(s0); CHKERRQ(ierr); 00708 PetscFunctionReturn(0); 00709 } 00710 00711 #undef __FUNCT__ 00712 #define __FUNCT__ "DMCreateFieldDecompositionDM_libMesh" 00713 static PetscErrorCode DMCreateFieldDecompositionDM_libMesh(DM dm, const char* ddesc, DM *ddm) 00714 { 00715 PetscFunctionBegin; 00716 PetscErrorCode ierr; 00717 PetscInt dtype, dcount, *dsizes; 00718 char ***dlists; 00719 PetscFunctionBegin; 00720 *ddm = PETSC_NULL; 00721 ierr = DMlibMeshParseDecompositionDescriptor_Private(dm,ddesc,&dtype,&dcount,&dsizes,&dlists); CHKERRQ(ierr); 00722 if(dtype == DMLIBMESH_FIELD_DECOMPOSITION){ 00723 ierr = DMlibMeshCreateFieldDecompositionDM(dm,dcount,dsizes,dlists,ddm); CHKERRQ(ierr); 00724 } 00725 else SETERRQ1(((PetscObject)dm)->comm, PETSC_ERR_PLIB, "Uexpected unknown decomposition type for field decomposition descriptor %s", ddesc); 00726 PetscFunctionReturn(0); 00727 } 00728 00729 #undef __FUNCT__ 00730 #define __FUNCT__ "DMCreateDomainDecompositionDM_libMesh" 00731 static PetscErrorCode DMCreateDomainDecompositionDM_libMesh(DM dm, const char* ddesc, DM *ddm) 00732 { 00733 PetscFunctionBegin; 00734 PetscErrorCode ierr; 00735 PetscInt dtype, dcount, *dsizes; 00736 char ***dlists; 00737 PetscFunctionBegin; 00738 *ddm = PETSC_NULL; 00739 ierr = DMlibMeshParseDecompositionDescriptor_Private(dm,ddesc,&dtype,&dcount,&dsizes,&dlists); CHKERRQ(ierr); 00740 if(dtype == DMLIBMESH_DOMAIN_DECOMPOSITION) { 00741 ierr = DMlibMeshCreateDomainDecompositionDM(dm,dcount,dsizes,dlists,ddm); CHKERRQ(ierr); 00742 } 00743 else SETERRQ1(((PetscObject)dm)->comm, PETSC_ERR_PLIB, "Uexpected unknown decomposition type for domain decomposition descriptor %s", ddesc); 00744 PetscFunctionReturn(0); 00745 } 00746 00747 #endif 00748 00749 00750 #undef __FUNCT__ 00751 #define __FUNCT__ "DMlibMeshFunction" 00752 static PetscErrorCode DMlibMeshFunction(DM dm, Vec x, Vec r) 00753 { 00754 PetscErrorCode ierr; 00755 PetscFunctionBegin; 00756 libmesh_assert(x); 00757 libmesh_assert(r); 00758 00759 NonlinearImplicitSystem* _sys; 00760 ierr = DMlibMeshGetSystem(dm, _sys);CHKERRQ(ierr); 00761 NonlinearImplicitSystem& sys = *_sys; 00762 PetscVector<Number>& X_sys = *libmesh_cast_ptr<PetscVector<Number>* >(sys.solution.get()); 00763 PetscVector<Number>& R_sys = *libmesh_cast_ptr<PetscVector<Number>* >(sys.rhs); 00764 PetscVector<Number> X_global(x, _sys->comm()), R(r, _sys->comm()); 00765 00766 // Use the systems update() to get a good local version of the parallel solution 00767 X_global.swap(X_sys); 00768 R.swap(R_sys); 00769 00770 _sys->get_dof_map().enforce_constraints_exactly(*_sys); 00771 _sys->update(); 00772 00773 // Swap back 00774 X_global.swap(X_sys); 00775 R.swap(R_sys); 00776 R.zero(); 00777 00778 // if the user has provided both function pointers and objects only the pointer 00779 // will be used, so catch that as an error 00780 if (_sys->nonlinear_solver->residual && _sys->nonlinear_solver->residual_object) 00781 libmesh_error_msg("ERROR: cannot specifiy both a function and object to compute the Residual!"); 00782 00783 if (_sys->nonlinear_solver->matvec && _sys->nonlinear_solver->residual_and_jacobian_object) 00784 libmesh_error_msg("ERROR: cannot specifiy both a function and object to compute the combined Residual & Jacobian!"); 00785 00786 if (_sys->nonlinear_solver->residual != NULL) 00787 _sys->nonlinear_solver->residual(*(_sys->current_local_solution.get()), R, *_sys); 00788 00789 else if (_sys->nonlinear_solver->residual_object != NULL) 00790 _sys->nonlinear_solver->residual_object->residual(*(_sys->current_local_solution.get()), R, *_sys); 00791 00792 else if (_sys->nonlinear_solver->matvec != NULL) 00793 _sys->nonlinear_solver->matvec(*(_sys->current_local_solution.get()), &R, NULL, *_sys); 00794 00795 else if (_sys->nonlinear_solver->residual_and_jacobian_object != NULL) 00796 _sys->nonlinear_solver->residual_and_jacobian_object->residual_and_jacobian(*(_sys->current_local_solution.get()), &R, NULL, *_sys); 00797 00798 else 00799 libmesh_error_msg("Error! Unable to compute residual and/or Jacobian!"); 00800 00801 R.close(); 00802 X_global.close(); 00803 PetscFunctionReturn(0); 00804 } 00805 00806 #if !PETSC_RELEASE_LESS_THAN(3,3,1) 00807 #undef __FUNCT__ 00808 #define __FUNCT__ "SNESFunction_DMlibMesh" 00809 static PetscErrorCode SNESFunction_DMlibMesh(SNES, Vec x, Vec r, void *ctx) 00810 { 00811 DM dm = (DM)ctx; 00812 PetscErrorCode ierr; 00813 PetscFunctionBegin; 00814 ierr = DMlibMeshFunction(dm,x,r);CHKERRQ(ierr); 00815 PetscFunctionReturn(0); 00816 } 00817 #endif 00818 00819 00820 #undef __FUNCT__ 00821 #define __FUNCT__ "DMlibMeshJacobian" 00822 #if PETSC_RELEASE_LESS_THAN(3,5,0) 00823 static PetscErrorCode DMlibMeshJacobian(DM dm, Vec x, Mat jac, Mat pc, MatStructure *msflag) 00824 #else 00825 static PetscErrorCode DMlibMeshJacobian(DM dm, Vec x, Mat jac, Mat pc) 00826 #endif 00827 { 00828 PetscErrorCode ierr; 00829 PetscFunctionBegin; 00830 NonlinearImplicitSystem* _sys; 00831 ierr = DMlibMeshGetSystem(dm, _sys); CHKERRQ(ierr); 00832 NonlinearImplicitSystem& sys = *_sys; 00833 00834 PetscMatrix<Number> the_pc(pc,sys.comm()); 00835 PetscMatrix<Number> Jac(jac,sys.comm()); 00836 PetscVector<Number>& X_sys = *libmesh_cast_ptr<PetscVector<Number>*>(sys.solution.get()); 00837 PetscMatrix<Number>& Jac_sys = *libmesh_cast_ptr<PetscMatrix<Number>*>(sys.matrix); 00838 PetscVector<Number> X_global(x, sys.comm()); 00839 00840 // Set the dof maps 00841 the_pc.attach_dof_map(sys.get_dof_map()); 00842 Jac.attach_dof_map(sys.get_dof_map()); 00843 00844 // Use the systems update() to get a good local version of the parallel solution 00845 X_global.swap(X_sys); 00846 Jac.swap(Jac_sys); 00847 00848 sys.get_dof_map().enforce_constraints_exactly(sys); 00849 sys.update(); 00850 00851 X_global.swap(X_sys); 00852 Jac.swap(Jac_sys); 00853 00854 the_pc.zero(); 00855 00856 // if the user has provided both function pointers and objects only the pointer 00857 // will be used, so catch that as an error 00858 if (sys.nonlinear_solver->jacobian && sys.nonlinear_solver->jacobian_object) 00859 libmesh_error_msg("ERROR: cannot specifiy both a function and object to compute the Jacobian!"); 00860 00861 if (sys.nonlinear_solver->matvec && sys.nonlinear_solver->residual_and_jacobian_object) 00862 libmesh_error_msg("ERROR: cannot specifiy both a function and object to compute the combined Residual & Jacobian!"); 00863 00864 if (sys.nonlinear_solver->jacobian != NULL) 00865 sys.nonlinear_solver->jacobian(*(sys.current_local_solution.get()), the_pc, sys); 00866 00867 else if (sys.nonlinear_solver->jacobian_object != NULL) 00868 sys.nonlinear_solver->jacobian_object->jacobian(*(sys.current_local_solution.get()), the_pc, sys); 00869 00870 else if (sys.nonlinear_solver->matvec != NULL) 00871 sys.nonlinear_solver->matvec(*(sys.current_local_solution.get()), NULL, &the_pc, sys); 00872 00873 else if (sys.nonlinear_solver->residual_and_jacobian_object != NULL) 00874 sys.nonlinear_solver->residual_and_jacobian_object->residual_and_jacobian(*(sys.current_local_solution.get()), NULL, &the_pc, sys); 00875 00876 else 00877 libmesh_error_msg("Error! Unable to compute residual and/or Jacobian!"); 00878 00879 the_pc.close(); 00880 Jac.close(); 00881 X_global.close(); 00882 #if PETSC_RELEASE_LESS_THAN(3,5,0) 00883 *msflag = SAME_NONZERO_PATTERN; 00884 #endif 00885 PetscFunctionReturn(0); 00886 } 00887 00888 #if !PETSC_RELEASE_LESS_THAN(3,3,1) 00889 #undef __FUNCT__ 00890 #define __FUNCT__ "SNESJacobian_DMlibMesh" 00891 #if PETSC_RELEASE_LESS_THAN(3,5,0) 00892 static PetscErrorCode SNESJacobian_DMlibMesh(SNES,Vec x,Mat *jac,Mat *pc, MatStructure* flag, void* ctx) 00893 #else 00894 static PetscErrorCode SNESJacobian_DMlibMesh(SNES,Vec x,Mat jac,Mat pc, void* ctx) 00895 #endif 00896 { 00897 DM dm = (DM)ctx; 00898 PetscErrorCode ierr; 00899 PetscFunctionBegin; 00900 #if PETSC_RELEASE_LESS_THAN(3,5,0) 00901 ierr = DMlibMeshJacobian(dm,x,*jac,*pc,flag); CHKERRQ(ierr); 00902 #else 00903 ierr = DMlibMeshJacobian(dm,x,jac,pc); CHKERRQ(ierr); 00904 #endif 00905 PetscFunctionReturn(0); 00906 } 00907 #endif 00908 00909 #undef __FUNCT__ 00910 #define __FUNCT__ "DMVariableBounds_libMesh" 00911 static PetscErrorCode DMVariableBounds_libMesh(DM dm, Vec xl, Vec xu) 00912 { 00913 PetscErrorCode ierr; 00914 NonlinearImplicitSystem* _sys; 00915 ierr = DMlibMeshGetSystem(dm, _sys); CHKERRQ(ierr); 00916 NonlinearImplicitSystem& sys = *_sys; 00917 PetscVector<Number> XL(xl, sys.comm()); 00918 PetscVector<Number> XU(xu, sys.comm()); 00919 PetscFunctionBegin; 00920 #if PETSC_VERSION_LESS_THAN(3,5,0) && PETSC_VERSION_RELEASE 00921 ierr = VecSet(xl, SNES_VI_NINF);CHKERRQ(ierr); 00922 ierr = VecSet(xu, SNES_VI_INF);CHKERRQ(ierr); 00923 #else 00924 ierr = VecSet(xl, PETSC_NINFINITY);CHKERRQ(ierr); 00925 ierr = VecSet(xu, PETSC_INFINITY);CHKERRQ(ierr); 00926 #endif 00927 if (sys.nonlinear_solver->bounds != NULL) 00928 sys.nonlinear_solver->bounds(XL,XU,sys); 00929 else if (sys.nonlinear_solver->bounds_object != NULL) 00930 sys.nonlinear_solver->bounds_object->bounds(XL,XU, sys); 00931 else 00932 SETERRQ(((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONG, "No bounds calculation in this libMesh object"); 00933 00934 PetscFunctionReturn(0); 00935 } 00936 00937 00938 #undef __FUNCT__ 00939 #define __FUNCT__ "DMCreateGlobalVector_libMesh" 00940 static PetscErrorCode DMCreateGlobalVector_libMesh(DM dm, Vec *x) 00941 { 00942 PetscFunctionBegin; 00943 PetscErrorCode ierr; 00944 DM_libMesh *dlm = (DM_libMesh *)(dm->data); 00945 PetscBool eq; 00946 00947 ierr = PetscObjectTypeCompare((PetscObject)dm, DMLIBMESH, &eq); CHKERRQ(ierr); 00948 00949 if (!eq) 00950 SETERRQ2(((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONG, "DM of type %s, not of type %s", ((PetscObject)dm)->type, DMLIBMESH); 00951 00952 if (!dlm->sys) 00953 SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONGSTATE, "No libMesh system set for DM_libMesh"); 00954 00955 NumericVector<Number>* nv = (dlm->sys->solution).get(); 00956 PetscVector<Number>* pv = dynamic_cast<PetscVector<Number>*>(nv); 00957 Vec v = pv->vec(); 00958 /* Unfortunately, currently this does not produce a ghosted vector, so nonlinear subproblem solves aren't going to be easily available. 00959 Should work fine for getting vectors out for linear subproblem solvers. */ 00960 if(dlm->embedding) { 00961 PetscInt n; 00962 ierr = VecCreate(((PetscObject)v)->comm, x); CHKERRQ(ierr); 00963 ierr = ISGetLocalSize(dlm->embedding, &n); CHKERRQ(ierr); 00964 ierr = VecSetSizes(*x,n,PETSC_DETERMINE); CHKERRQ(ierr); 00965 ierr = VecSetType(*x,((PetscObject)v)->type_name); CHKERRQ(ierr); 00966 ierr = VecSetFromOptions(*x); CHKERRQ(ierr); 00967 ierr = VecSetUp(*x); CHKERRQ(ierr); 00968 } 00969 else { 00970 ierr = VecDuplicate(v,x); CHKERRQ(ierr); 00971 } 00972 ierr = PetscObjectCompose((PetscObject)*x,"DM",(PetscObject)dm); CHKERRQ(ierr); 00973 PetscFunctionReturn(0); 00974 } 00975 00976 00977 00978 00979 #undef __FUNCT__ 00980 #define __FUNCT__ "DMCreateMatrix_libMesh" 00981 #if PETSC_VERSION_LT(3,5,0) 00982 static PetscErrorCode DMCreateMatrix_libMesh(DM dm, const MatType, Mat *A) 00983 #else 00984 static PetscErrorCode DMCreateMatrix_libMesh(DM dm, Mat *A) 00985 #endif 00986 { 00987 PetscFunctionBegin; 00988 PetscErrorCode ierr; 00989 DM_libMesh *dlm = (DM_libMesh *)(dm->data); 00990 PetscBool eq; 00991 00992 ierr = PetscObjectTypeCompare((PetscObject)dm, DMLIBMESH, &eq); CHKERRQ(ierr); 00993 00994 if (!eq) 00995 SETERRQ2(((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONG, "DM of type %s, not of type %s", ((PetscObject)dm)->type, DMLIBMESH); 00996 00997 if (!dlm->sys) 00998 SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONGSTATE, "No libMesh system set for DM_libMesh"); 00999 01000 *A = (dynamic_cast<PetscMatrix<Number>*>(dlm->sys->matrix))->mat(); 01001 ierr = PetscObjectReference((PetscObject)(*A)); CHKERRQ(ierr); 01002 PetscFunctionReturn(0); 01003 } 01004 01005 01006 #undef __FUNCT__ 01007 #define __FUNCT__ "DMView_libMesh" 01008 static PetscErrorCode DMView_libMesh(DM dm, PetscViewer viewer) 01009 { 01010 PetscErrorCode ierr; 01011 PetscBool isascii; 01012 const char *name, *prefix; 01013 DM_libMesh *dlm = (DM_libMesh*)dm->data; 01014 PetscFunctionBegin; 01015 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii); CHKERRQ(ierr); 01016 if(isascii) { 01017 ierr = PetscObjectGetName((PetscObject)dm, &name); CHKERRQ(ierr); 01018 ierr = PetscObjectGetOptionsPrefix((PetscObject)dm, &prefix); CHKERRQ(ierr); 01019 ierr = PetscViewerASCIIPrintf(viewer, "DM libMesh with name %s and prefix %s\n", name, prefix); CHKERRQ(ierr); 01020 ierr = PetscViewerASCIIPrintf(viewer, "blocks:", name, prefix); CHKERRQ(ierr); 01021 std::map<std::string,unsigned int>::iterator bit = dlm->blockids->begin(); 01022 std::map<std::string,unsigned int>::const_iterator bend = dlm->blockids->end(); 01023 for(; bit != bend; ++bit) { 01024 ierr = PetscViewerASCIIPrintf(viewer, "(%s,%D) ", bit->first.c_str(), bit->second); CHKERRQ(ierr); 01025 } 01026 ierr = PetscViewerASCIIPrintf(viewer, "\n"); CHKERRQ(ierr); 01027 ierr = PetscViewerASCIIPrintf(viewer, "variables:", name, prefix); CHKERRQ(ierr); 01028 std::map<std::string,unsigned int>::iterator vit = dlm->varids->begin(); 01029 std::map<std::string,unsigned int>::const_iterator vend = dlm->varids->end(); 01030 for(; vit != vend; ++vit) { 01031 ierr = PetscViewerASCIIPrintf(viewer, "(%s,%D) ", vit->first.c_str(), vit->second); CHKERRQ(ierr); 01032 } 01033 ierr = PetscViewerASCIIPrintf(viewer, "\n"); CHKERRQ(ierr); 01034 if(dlm->decomposition_type == DMLIBMESH_NO_DECOMPOSITION) { 01035 ierr = PetscViewerASCIIPrintf(viewer, "No decomposition\n"); CHKERRQ(ierr); 01036 } 01037 else { 01038 if(dlm->decomposition_type == DMLIBMESH_FIELD_DECOMPOSITION) { 01039 ierr = PetscViewerASCIIPrintf(viewer, "Field decomposition by variable: "); CHKERRQ(ierr); 01040 } 01041 else if(dlm->decomposition_type == DMLIBMESH_DOMAIN_DECOMPOSITION) { 01042 ierr = PetscViewerASCIIPrintf(viewer, "Domain decomposition by block: "); CHKERRQ(ierr); 01043 } 01044 else SETERRQ1(((PetscObject)dm)->comm, PETSC_ERR_PLIB, "Unexpected decomposition type: %D", dlm->decomposition_type); 01045 /* FIX: decompositions might have different sizes and components on different ranks. */ 01046 for(unsigned int d = 0; d < dlm->decomposition->size(); ++d) { 01047 std::set<unsigned int>::iterator dbegin = (*dlm->decomposition)[d].begin(); 01048 std::set<unsigned int>::iterator dit = (*dlm->decomposition)[d].begin(); 01049 std::set<unsigned int>::iterator dend = (*dlm->decomposition)[d].end(); 01050 for(; dit != dend; ++dit) { 01051 if(dit != dbegin) { 01052 ierr = PetscViewerASCIIPrintf(viewer, ","); CHKERRQ(ierr); 01053 } 01054 ierr = PetscViewerASCIIPrintf(viewer, "%D", *dit); CHKERRQ(ierr); 01055 } 01056 ierr = PetscViewerASCIIPrintf(viewer, ";"); CHKERRQ(ierr); 01057 } 01058 ierr = PetscViewerASCIIPrintf(viewer, "\n"); CHKERRQ(ierr); 01059 } 01060 } 01061 01062 PetscFunctionReturn(0); 01063 } 01064 01065 #undef __FUNCT__ 01066 #define __FUNCT__ "DMSetUp_libMesh" 01067 static PetscErrorCode DMSetUp_libMesh(DM dm) 01068 { 01069 PetscFunctionBegin; 01070 PetscErrorCode ierr; 01071 DM_libMesh *dlm = (DM_libMesh *)(dm->data); 01072 PetscBool eq; 01073 01074 ierr = PetscObjectTypeCompare((PetscObject)dm, DMLIBMESH, &eq); CHKERRQ(ierr); 01075 01076 if (!eq) 01077 SETERRQ2(((PetscObject)dm)->comm, PETSC_ERR_ARG_WRONG, "DM of type %s, not of type %s", ((PetscObject)dm)->type, DMLIBMESH); 01078 01079 if (!dlm->sys) 01080 SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONGSTATE, "No libMesh system set for DM_libMesh"); 01081 /* 01082 Do not evaluate function, Jacobian or bounds for an embedded DM -- the subproblem might not have enough information for that. 01083 */ 01084 if(!dlm->embedding) { 01085 #if PETSC_RELEASE_LESS_THAN(3,3,1) 01086 ierr = DMSetFunction(dm, DMlibMeshFunction); CHKERRQ(ierr); 01087 ierr = DMSetJacobian(dm, DMlibMeshJacobian); CHKERRQ(ierr); 01088 #else 01089 ierr = DMSNESSetFunction(dm, SNESFunction_DMlibMesh, (void*)dm); CHKERRQ(ierr); 01090 ierr = DMSNESSetJacobian(dm, SNESJacobian_DMlibMesh, (void*)dm); CHKERRQ(ierr); 01091 #endif 01092 if (dlm->sys->nonlinear_solver->bounds || dlm->sys->nonlinear_solver->bounds_object) 01093 ierr = DMSetVariableBounds(dm, DMVariableBounds_libMesh); CHKERRQ(ierr); 01094 } 01095 else { 01096 /* 01097 Fow now we don't implement even these, although a linear "Dirichlet" subproblem is well-defined. 01098 Creating the submatrix, however, might require extracting the submatrix preallocation from an unassembled matrix. 01099 */ 01100 dm->ops->createglobalvector = 0; 01101 dm->ops->creatematrix = 0; 01102 } 01103 PetscFunctionReturn(0); 01104 } 01105 01106 01107 01108 01109 #undef __FUNCT__ 01110 #define __FUNCT__ "DMDestroy_libMesh" 01111 static PetscErrorCode DMDestroy_libMesh(DM dm) 01112 { 01113 DM_libMesh *dlm = (DM_libMesh*)(dm->data); 01114 PetscErrorCode ierr; 01115 PetscFunctionBegin; 01116 delete dlm->varids; 01117 delete dlm->varnames; 01118 delete dlm->blockids; 01119 delete dlm->blocknames; 01120 delete dlm->decomposition; 01121 ierr = ISDestroy(&dlm->embedding); CHKERRQ(ierr); 01122 ierr = PetscFree(dm->data); CHKERRQ(ierr); 01123 01124 PetscFunctionReturn(0); 01125 } 01126 01127 EXTERN_C_BEGIN 01128 #undef __FUNCT__ 01129 #define __FUNCT__ "DMCreate_libMesh" 01130 PetscErrorCode DMCreate_libMesh(DM dm) 01131 { 01132 PetscErrorCode ierr; 01133 DM_libMesh *dlm; 01134 01135 PetscFunctionBegin; 01136 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 01137 #if PETSC_RELEASE_LESS_THAN(3,5,0) 01138 ierr = PetscNewLog(dm,DM_libMesh,&dlm);CHKERRQ(ierr); 01139 #else 01140 ierr = PetscNewLog(dm,&dlm);CHKERRQ(ierr); 01141 #endif 01142 dm->data = dlm; 01143 01144 /* DMlibMesh impl */ 01145 dlm->varids = new(std::map<std::string, unsigned int>); 01146 dlm->blockids = new(std::map<std::string, unsigned int>); 01147 dlm->varnames = new(std::map<unsigned int, std::string>); 01148 dlm->blocknames = new(std::map<unsigned int, std::string>); 01149 dlm->decomposition = PETSC_NULL; 01150 dlm->decomposition_type = DMLIBMESH_NO_DECOMPOSITION; 01151 01152 /* DM API */ 01153 dm->ops->createglobalvector = DMCreateGlobalVector_libMesh; 01154 dm->ops->createlocalvector = 0; // DMCreateLocalVector_libMesh; 01155 dm->ops->getcoloring = 0; // DMGetColoring_libMesh; 01156 dm->ops->creatematrix = DMCreateMatrix_libMesh; 01157 dm->ops->createinterpolation= 0; // DMCreateInterpolation_libMesh; 01158 01159 dm->ops->refine = 0; // DMRefine_libMesh; 01160 dm->ops->coarsen = 0; // DMCoarsen_libMesh; 01161 dm->ops->getinjection = 0; // DMGetInjection_libMesh; 01162 dm->ops->getaggregates = 0; // DMGetAggregates_libMesh; 01163 01164 #if PETSC_RELEASE_LESS_THAN(3,3,1) 01165 dm->ops->createfielddecompositiondm = DMCreateFieldDecompositionDM_libMesh; 01166 dm->ops->createdomaindecompositiondm = DMCreateDomainDecompositionDM_libMesh; 01167 #endif 01168 dm->ops->createfielddecomposition = DMCreateFieldDecomposition_libMesh; 01169 dm->ops->createdomaindecomposition = DMCreateDomainDecomposition_libMesh; 01170 01171 dm->ops->destroy = DMDestroy_libMesh; 01172 dm->ops->view = DMView_libMesh; 01173 dm->ops->setfromoptions = 0; // DMSetFromOptions_libMesh; 01174 dm->ops->setup = DMSetUp_libMesh; 01175 01176 /* DMlibMesh API */ 01177 #if PETSC_RELEASE_LESS_THAN(3,4,0) 01178 ierr = PetscObjectComposeFunction((PetscObject)dm,"DMlibMeshSetSystem_C",PETSC_NULL,(PetscVoidFunction)DMlibMeshSetSystem_libMesh);CHKERRQ(ierr); 01179 ierr = PetscObjectComposeFunction((PetscObject)dm,"DMlibMeshGetSystem_C",PETSC_NULL,(PetscVoidFunction)DMlibMeshGetSystem_libMesh);CHKERRQ(ierr); 01180 #else 01181 ierr = PetscObjectComposeFunction((PetscObject)dm,"DMlibMeshSetSystem_C",DMlibMeshSetSystem_libMesh);CHKERRQ(ierr); 01182 ierr = PetscObjectComposeFunction((PetscObject)dm,"DMlibMeshGetSystem_C",DMlibMeshGetSystem_libMesh);CHKERRQ(ierr); 01183 #endif 01184 01185 PetscFunctionReturn(0); 01186 } 01187 EXTERN_C_END 01188 01189 01190 #endif // #if !PETSC_VERSION_LESS_THAN(3,3,0)