]> SALOME platform Git repositories - tools/solverlab.git/commitdiff
Salome HOME
Tested mesh fast equivalence when giving an input field + MPI parallelisation of...
authormichael <michael@localhost.localdomain>
Fri, 10 Dec 2021 13:34:48 +0000 (14:34 +0100)
committermichael <michael@localhost.localdomain>
Fri, 10 Dec 2021 13:34:48 +0000 (14:34 +0100)
CoreFlows/Models/inc/DiffusionEquation.hxx
CoreFlows/Models/src/DiffusionEquation.cxx

index ec2692c3d4207dc13fa152afcbbac21284be77c1..16e24c505b32671a404c33a70be314798cf188d8 100755 (executable)
@@ -106,17 +106,9 @@ public :
        vector<string> getInputFieldsNames();
        void setInputField(const string& nameField, Field& inputField );//supply of a required input field
 
-       void setFluidTemperatureField(Field coupledTemperatureField){
-               coupledTemperatureField.getMesh().checkFastEquivalWith(_mesh);
-               _fluidTemperatureField=coupledTemperatureField;
-               _fluidTemperatureFieldSet=true;
-       };
-       void setFluidTemperature(double fluidTemperature){
-       _fluidTemperature=fluidTemperature;
-       }
-       Field& getFluidTemperatureField(){
-               return _fluidTemperatureField;
-       }
+       void setFluidTemperatureField(Field coupledTemperatureField);
+       void setFluidTemperature(double fluidTemperature){      _fluidTemperature=fluidTemperature;     }
+       Field& getFluidTemperatureField(){  return _fluidTemperatureField;      }
        
        /*** get output fields names for postprocessing or coupling ***/
        vector<string> getOutputFieldsNames() ;//liste tous les champs que peut fournir le code pour le postraitement
index 4c4a7181251067a8af6c5bb6bcb9e21f6997d227..f76987032792459d32ef290c8241482c1f586d28 100755 (executable)
@@ -269,10 +269,35 @@ void DiffusionEquation::initialize()
        //Linear solver
        KSPCreate(PETSC_COMM_WORLD, &_ksp);
        KSPSetType(_ksp, _ksptype);
-       // if(_ksptype == KSPGMRES) KSPGMRESSetRestart(_ksp,10000);
-       KSPSetTolerances(_ksp,_precision,_precision,PETSC_DEFAULT,_maxPetscIts);
        KSPGetPC(_ksp, &_pc);
-       PCSetType(_pc, _pctype);
+       if(_mpi_size==1 )
+               PCSetType(_pc, _pctype);
+       else
+       {
+               PCSetType(_pc, PCBJACOBI);//Global preconditioner is block jacobi
+               if(_pctype != (char*)&PCILU)//Default pc type is ilu
+               {
+                       PetscOptionsSetValue(NULL,"-sub_pc_type ",_pctype);
+                       PetscOptionsSetValue(NULL,"-sub_ksp_type ","preonly");
+                       //If the above setvalue does not work, try the following
+                       /*
+                       KSPSetUp(_ksp);//to set the block Jacobi data structures (including creation of an internal KSP context for each block)
+                       KSP * subKSP;
+                       PC subpc;
+                       int nlocal;//nb local blocs (should equal 1)
+                       PCBJacobiGetSubKSP(_pc,&nlocal,NULL,&subKSP);
+                       if(nlocal==1)
+                       {
+                               KSPSetType(subKSP[0], KSPPREONLY);//local block solver is same as global
+                               KSPGetPC(subKSP[0],&subpc);
+                               PCSetType(subpc,_pctype);
+                       }
+                       else
+                               throw CdmathException("PC Block Jacobi, more than one block in this processor!!");
+                       */ 
+               }
+       }
+       KSPSetTolerances(_ksp,_precision,_precision,PETSC_DEFAULT,_maxPetscIts);
 
        _initializedMemory=true;
        save();//save initial data
@@ -890,6 +915,9 @@ DiffusionEquation::getOutputField(const string& nameField )
 void
 DiffusionEquation::setInputField(const string& nameField, Field& inputField )
 {
+       if(!_initialDataSet)
+               throw CdmathException("!!!!!!!! DiffusionEquation::setInputField set initial field first");
+
        if(nameField=="FluidTemperature" || nameField=="FLUIDTEMPERATURE" || nameField=="TemperatureFluide" || nameField=="TEMPERATUREFLUIDE")
                return setFluidTemperatureField( inputField) ;
        else if(nameField=="HeatPower" || nameField=="HEATPOWER" || nameField=="PuissanceThermique" || nameField=="PUISSANCETHERMIQUE" )
@@ -900,3 +928,13 @@ DiffusionEquation::setInputField(const string& nameField, Field& inputField )
         throw CdmathException("DiffusionEquation::setInputField error : Unknown Field name");
     }
 }
+
+void 
+DiffusionEquation::setFluidTemperatureField(Field coupledTemperatureField){
+       if(!_initialDataSet)
+               throw CdmathException("!!!!!!!! DiffusionEquation::setFluidTemperatureField set initial field first");
+
+       coupledTemperatureField.getMesh().checkFastEquivalWith(_mesh);
+       _fluidTemperatureField=coupledTemperatureField;
+       _fluidTemperatureFieldSet=true;
+};