Back to home page

MITgcm

 
 

    


File indexing completed on 2023-05-28 05:10:57 UTC

view on githubraw file Latest commit b4daa243 on 2023-05-28 03:53:22 UTC
b4daa24319 Shre*0001 /*
                0002 ##########################################################
                0003 # This file is part of the AdjoinableMPI library         #
                0004 # released under the MIT License.                        #
                0005 # The full COPYRIGHT notice can be found in the top      #
                0006 # level directory of the AdjoinableMPI distribution.     #
                0007 ########################################################## 
                0008 */
                0009 
                0010 /* The following to set the AMPI_FORTRANCOMPATIBLE
                0011  * if it was placed in the configure command. */
                0012 #include "ampi/userIF/libConfig.h"
                0013 
                0014 #include <stdlib.h>
                0015 #include <stdio.h>
                0016 #include <string.h>
                0017 #include <assert.h>
                0018 #include "ampi/adTool/support.h"
                0019 #include <mpi.h>
                0020 #include "ampi/ampi.h"
                0021 #include "ampi/libCommon/modified.h"
                0022 
                0023 struct AMPI_ShadowComm_list {
                0024   struct AMPI_ShadowComm_list *next_p;
                0025   MPI_Comm comm ;
                0026   MPI_Comm shadowComm ;
                0027 } ;
                0028 
                0029 struct AMPI_ShadowComm_list * ADTOOL_AMPI_SHADOWCOMM_LIST = NULL ;
                0030 
                0031 struct AMPI_Request_stack {
                0032   struct AMPI_Request_stack *next_p;
                0033   void *buf ;
                0034   void *adjointBuf ;
                0035   int count ;
                0036   MPI_Datatype datatype ;
                0037   int endPoint ;
                0038   int tag ;
                0039   enum AMPI_PairedWith_E pairedWith;
                0040   MPI_Comm comm;
                0041   enum AMPI_Activity_E isActive;
                0042   enum AMPI_Request_origin_E origin;
                0043 } ;
                0044 
                0045 int AMPI_Init_NT(int* argc, char*** argv) {
                0046   int rc = MPI_Init(argc, argv);
                0047   ADTOOL_AMPI_setupTypes() ;
                0048   MPI_Comm worldDup ;
                0049   int rc2 = MPI_Comm_dup(MPI_COMM_WORLD, &worldDup) ;
                0050   assert(rc2==MPI_SUCCESS);
                0051   ADTOOL_AMPI_SHADOWCOMM_LIST = NULL ;
                0052   ADTOOL_AMPI_addShadowComm(MPI_COMM_WORLD, worldDup) ;
                0053   ourADTOOL_AMPI_FPCollection.pushBcastInfo_fp=&ADTOOL_AMPI_pushBcastInfo;
                0054   ourADTOOL_AMPI_FPCollection.popBcastInfo_fp=&ADTOOL_AMPI_popBcastInfo;
                0055   ourADTOOL_AMPI_FPCollection.pushDoubleArray_fp=&ADTOOL_AMPI_pushDoubleArray;
                0056   ourADTOOL_AMPI_FPCollection.popDoubleArray_fp=&ADTOOL_AMPI_popDoubleArray;
                0057   ourADTOOL_AMPI_FPCollection.pushReduceInfo_fp=&ADTOOL_AMPI_pushReduceInfo; 
                0058   ourADTOOL_AMPI_FPCollection.popReduceCountAndType_fp=&ADTOOL_AMPI_popReduceCountAndType;
                0059   ourADTOOL_AMPI_FPCollection.popReduceInfo_fp=&ADTOOL_AMPI_popReduceInfo; 
                0060   ourADTOOL_AMPI_FPCollection.pushSRinfo_fp=&ADTOOL_AMPI_pushSRinfo;
                0061   ourADTOOL_AMPI_FPCollection.popSRinfo_fp=&ADTOOL_AMPI_popSRinfo;
                0062   ourADTOOL_AMPI_FPCollection.pushGSinfo_fp=&ADTOOL_AMPI_pushGSinfo;
                0063   ourADTOOL_AMPI_FPCollection.popGScommSizeForRootOrNull_fp=&ADTOOL_AMPI_popGScommSizeForRootOrNull;
                0064   ourADTOOL_AMPI_FPCollection.popGSinfo_fp=&ADTOOL_AMPI_popGSinfo;
                0065   ourADTOOL_AMPI_FPCollection.pushGSVinfo_fp=&ADTOOL_AMPI_pushGSVinfo;
                0066   ourADTOOL_AMPI_FPCollection.popGSVinfo_fp=&ADTOOL_AMPI_popGSVinfo;
                0067   ourADTOOL_AMPI_FPCollection.push_CallCode_fp=&ADTOOL_AMPI_push_CallCode;
                0068   ourADTOOL_AMPI_FPCollection.pop_CallCode_fp=&ADTOOL_AMPI_pop_CallCode;
                0069   ourADTOOL_AMPI_FPCollection.push_AMPI_Request_fp=&ADTOOL_AMPI_push_AMPI_Request;
                0070   ourADTOOL_AMPI_FPCollection.pop_AMPI_Request_fp=&ADTOOL_AMPI_pop_AMPI_Request;
                0071   ourADTOOL_AMPI_FPCollection.push_request_fp=&ADTOOL_AMPI_push_request;
                0072   ourADTOOL_AMPI_FPCollection.pop_request_fp=&ADTOOL_AMPI_pop_request;
                0073   ourADTOOL_AMPI_FPCollection.push_comm_fp=&ADTOOL_AMPI_push_comm;
                0074   ourADTOOL_AMPI_FPCollection.pop_comm_fp=&ADTOOL_AMPI_pop_comm;
                0075   ourADTOOL_AMPI_FPCollection.rawData_fp=&ADTOOL_AMPI_rawData;
                0076   ourADTOOL_AMPI_FPCollection.rawDataV_fp=&ADTOOL_AMPI_rawDataV;
                0077   ourADTOOL_AMPI_FPCollection.packDType_fp=&ADTOOL_AMPI_packDType;
                0078   ourADTOOL_AMPI_FPCollection.unpackDType_fp=&ADTOOL_AMPI_unpackDType;
                0079   ourADTOOL_AMPI_FPCollection.writeData_fp=&ADTOOL_AMPI_writeData;
                0080   ourADTOOL_AMPI_FPCollection.writeDataV_fp=&ADTOOL_AMPI_writeDataV;
                0081   ourADTOOL_AMPI_FPCollection.rawAdjointData_fp=&ADTOOL_AMPI_rawAdjointData;
                0082   ourADTOOL_AMPI_FPCollection.Turn_fp=&ADTOOL_AMPI_Turn;
                0083   ourADTOOL_AMPI_FPCollection.mapBufForAdjoint_fp=&ADTOOL_AMPI_mapBufForAdjoint;
                0084   ourADTOOL_AMPI_FPCollection.setBufForAdjoint_fp=&ADTOOL_AMPI_setBufForAdjoint;
                0085   ourADTOOL_AMPI_FPCollection.getAdjointCount_fp=&ADTOOL_AMPI_getAdjointCount;
                0086   ourADTOOL_AMPI_FPCollection.setAdjointCount_fp=&ADTOOL_AMPI_setAdjointCount;
                0087   ourADTOOL_AMPI_FPCollection.setAdjointCountAndTempBuf_fp=&ADTOOL_AMPI_setAdjointCountAndTempBuf;
                0088   ourADTOOL_AMPI_FPCollection.allocateTempBuf_fp=&ADTOOL_AMPI_allocateTempBuf;
                0089   ourADTOOL_AMPI_FPCollection.releaseAdjointTempBuf_fp=&ADTOOL_AMPI_releaseAdjointTempBuf;
                0090   ourADTOOL_AMPI_FPCollection.adjointMultiply_fp=&ADTOOL_AMPI_adjointMultiply;
                0091   ourADTOOL_AMPI_FPCollection.adjointMin_fp=&ADTOOL_AMPI_adjointMin;
                0092   ourADTOOL_AMPI_FPCollection.adjointMax_fp=&ADTOOL_AMPI_adjointMax;
                0093   ourADTOOL_AMPI_FPCollection.multiplyAdjoint_fp=&ADTOOL_AMPI_multiplyAdjoint;
                0094   ourADTOOL_AMPI_FPCollection.divideAdjoint_fp=&ADTOOL_AMPI_divideAdjoint;
                0095   ourADTOOL_AMPI_FPCollection.equalAdjoints_fp=&ADTOOL_AMPI_equalAdjoints;
                0096   ourADTOOL_AMPI_FPCollection.incrementAdjoint_fp=&ADTOOL_AMPI_incrementAdjoint;
                0097   ourADTOOL_AMPI_FPCollection.nullifyAdjoint_fp=&ADTOOL_AMPI_nullifyAdjoint;
                0098   ourADTOOL_AMPI_FPCollection.setupTypes_fp=&ADTOOL_AMPI_setupTypes;
                0099   ourADTOOL_AMPI_FPCollection.cleanupTypes_fp=&ADTOOL_AMPI_cleanupTypes;
                0100   ourADTOOL_AMPI_FPCollection.FW_rawType_fp=&ADTOOL_AMPI_FW_rawType;
                0101   ourADTOOL_AMPI_FPCollection.BW_rawType_fp=&ADTOOL_AMPI_BW_rawType;
                0102   ourADTOOL_AMPI_FPCollection.isActiveType_fp=&ADTOOL_AMPI_isActiveType;
                0103   ourADTOOL_AMPI_FPCollection.allocateTempActiveBuf_fp=&ADTOOL_AMPI_allocateTempActiveBuf;
                0104   ourADTOOL_AMPI_FPCollection.releaseTempActiveBuf_fp=&ADTOOL_AMPI_releaseTempActiveBuf;
                0105   ourADTOOL_AMPI_FPCollection.copyActiveBuf_fp=&ADTOOL_AMPI_copyActiveBuf;
                0106   ourADTOOL_AMPI_FPCollection.tangentMultiply_fp=&ADTOOL_AMPI_tangentMultiply ;
                0107   ourADTOOL_AMPI_FPCollection.tangentMin_fp=&ADTOOL_AMPI_tangentMin ;
                0108   ourADTOOL_AMPI_FPCollection.tangentMax_fp=&ADTOOL_AMPI_tangentMax ;
                0109   ourADTOOL_AMPI_FPCollection.pushBuffer_fp=&ADTOOL_AMPI_pushBuffer ;
                0110   ourADTOOL_AMPI_FPCollection.popBuffer_fp=&ADTOOL_AMPI_popBuffer ;
                0111   ourADTOOL_AMPI_FPCollection.addShadowComm_fp=&ADTOOL_AMPI_addShadowComm ;
                0112   ourADTOOL_AMPI_FPCollection.getShadowComm_fp=&ADTOOL_AMPI_getShadowComm ;
                0113   ourADTOOL_AMPI_FPCollection.delShadowComm_fp=&ADTOOL_AMPI_delShadowComm ;
                0114 #ifdef AMPI_FORTRANCOMPATIBLE
                0115   ourADTOOL_AMPI_FPCollection.fortransetuptypes__fp=&adtool_ampi_fortransetuptypes_;
                0116   ourADTOOL_AMPI_FPCollection.fortrancleanuptypes__fp=&adtool_ampi_fortrancleanuptypes_;
                0117 #endif
                0118   return rc ;
                0119 }
                0120 
                0121 static struct AMPI_Request_stack* requestStackTop=0 ;
                0122 void ADTOOL_AMPI_pushBcastInfo(void* buf,
                0123                    int count,
                0124                    MPI_Datatype datatype,
                0125                    int root,
                0126                    MPI_Comm comm) {
                0127 }
                0128 
                0129 void ADTOOL_AMPI_popBcastInfo(void** buf,
                0130                   int* count,
                0131                   MPI_Datatype* datatype,
                0132                   int* root,
                0133                   MPI_Comm* comm,
                0134                   void **idx) {
                0135 }
                0136 
                0137 void ADTOOL_AMPI_pushDoubleArray(void* buf,
                0138                  int count) {
                0139 }
                0140 
                0141 void ADTOOL_AMPI_popDoubleArray(double* buf,
                0142                 int* count) {
                0143 }
                0144 
                0145 void ADTOOL_AMPI_pushReduceInfo(void* sbuf,
                0146                 void* rbuf,
                0147                 void* resultData,
                0148                 int pushResultData, /* push resultData if true */
                0149                 int count,
                0150                 MPI_Datatype datatype,
                0151                 MPI_Op op,
                0152                 int root,
                0153                 MPI_Comm comm) {
                0154 }
                0155 
                0156 void ADTOOL_AMPI_popReduceCountAndType(int* count,
                0157                        MPI_Datatype* datatype) {
                0158 }
                0159 
                0160 void ADTOOL_AMPI_popReduceInfo(void** sbuf,
                0161                    void** rbuf,
                0162                    void** prevData,
                0163                    void** resultData,
                0164                    int* count,
                0165                    MPI_Op* op,
                0166                    int* root,
                0167                    MPI_Comm* comm,
                0168                    void **idx) {
                0169 }
                0170 
                0171 extern void pushInteger4Array(int *x, int n) ;
                0172 extern void popInteger4Array(int *x, int n) ;
                0173 
                0174 void ADTOOL_AMPI_pushSRinfo(void* buf, 
                0175                 int count,
                0176                 MPI_Datatype datatype, 
                0177                 int src, 
                0178                 int tag,
                0179                 AMPI_PairedWith pairedWith,
                0180                 MPI_Comm comm) {
                0181   /* [llh] TODO: this is not nice: we should call pushInteger4() instead ! */
                0182   pushInteger4Array(&src,1) ;
                0183   pushInteger4Array(&tag,1) ;
                0184 }
                0185 
                0186 void ADTOOL_AMPI_popSRinfo(void** buf,
                0187                int* count,
                0188                MPI_Datatype* datatype, 
                0189                int* src, 
                0190                int* tag,
                0191                AMPI_PairedWith* pairedWith,
                0192                MPI_Comm* comm,
                0193                void **idx) { 
                0194   /* [llh] TODO: this is not nice: we should call popInteger4() instead ! */
                0195   popInteger4Array(tag,1) ;
                0196   popInteger4Array(src,1) ;
                0197 }
                0198 
                0199 void ADTOOL_AMPI_pushGSinfo(int commSizeForRootOrNull,
                0200                             void *rbuf,
                0201                             int rcnt,
                0202                             MPI_Datatype rtype,
                0203                             void *buf,
                0204                             int count,
                0205                             MPI_Datatype type,
                0206                             int  root,
                0207                             MPI_Comm comm) {
                0208   pushInteger4Array(&commSizeForRootOrNull,1) ;
                0209 }
                0210 
                0211 void ADTOOL_AMPI_popGScommSizeForRootOrNull(int *commSizeForRootOrNull) {
                0212   popInteger4Array(commSizeForRootOrNull,1) ;
                0213 }
                0214 
                0215 void ADTOOL_AMPI_popGSinfo(int commSizeForRootOrNull,
                0216                            void **rbuf,
                0217                            int *rcnt,
                0218                            MPI_Datatype *rtype,
                0219                            void **buf,
                0220                            int *count,
                0221                            MPI_Datatype *type,
                0222                            int *root,
                0223                            MPI_Comm *comm) {
                0224 }
                0225 
                0226 void ADTOOL_AMPI_pushGSVinfo(int commSizeForRootOrNull,
                0227                              void *rbuf,
                0228                              int *rcnts,
                0229                              int *displs,
                0230                              MPI_Datatype rtype,
                0231                              void *buf,
                0232                              int  count,
                0233                              MPI_Datatype type,
                0234                              int  root,
                0235                              MPI_Comm comm) {
                0236   pushInteger4Array(&commSizeForRootOrNull,1) ;
                0237 }
                0238 
                0239 void ADTOOL_AMPI_popGSVinfo(int commSizeForRootOrNull,
                0240                             void **rbuf,
                0241                             int *rcnts,
                0242                             int *displs,
                0243                             MPI_Datatype *rtype,
                0244                             void **buf,
                0245                             int *count,
                0246                             MPI_Datatype *type,
                0247                             int *root,
                0248                             MPI_Comm *comm) {
                0249 }
                0250 
                0251 void ADTOOL_AMPI_push_CallCode(enum AMPI_CallCode_E thisCall) { 
                0252 }
                0253 
                0254 void ADTOOL_AMPI_pop_CallCode(enum AMPI_CallCode_E *thisCall) { 
                0255 }
                0256 
                0257 void ADTOOL_AMPI_push_AMPI_Request(struct AMPI_Request_S  *ampiRequest) { 
                0258   struct AMPI_Request_stack* newTop =
                0259     (struct AMPI_Request_stack*)malloc(sizeof(struct AMPI_Request_stack)) ;
                0260   newTop->next_p = requestStackTop ;
                0261   newTop->buf = ampiRequest->buf ;
                0262   newTop->adjointBuf = ampiRequest->adjointBuf ;
                0263   newTop->count = ampiRequest->count ;
                0264   newTop->datatype = ampiRequest->datatype ;
                0265   newTop->endPoint = ampiRequest->endPoint ;
                0266   newTop->tag = ampiRequest->tag ;
                0267   newTop->pairedWith = ampiRequest->pairedWith ;
                0268   newTop->comm = ampiRequest->comm ;
                0269   newTop->origin = ampiRequest->origin ;
                0270   requestStackTop = newTop ;
                0271 }
                0272 
                0273 void ADTOOL_AMPI_pop_AMPI_Request(struct AMPI_Request_S  *ampiRequest) { 
                0274   struct AMPI_Request_stack* oldTop = requestStackTop ;
                0275   ampiRequest->buf = oldTop->buf ;
                0276   ampiRequest->adjointBuf = oldTop->adjointBuf ;
                0277   ampiRequest->count = oldTop->count ;
                0278   ampiRequest->datatype = oldTop->datatype ;
                0279   ampiRequest->endPoint = oldTop->endPoint ;
                0280   ampiRequest->tag = oldTop->tag ;
                0281   ampiRequest->pairedWith = oldTop->pairedWith ;
                0282   ampiRequest->comm = oldTop->comm ;
                0283   ampiRequest->origin = oldTop->origin ;
                0284   requestStackTop = oldTop->next_p ;
                0285   free(oldTop) ;
                0286 }
                0287 
                0288 void ADTOOL_AMPI_push_request(MPI_Request request) { 
                0289 } 
                0290 
                0291 MPI_Request ADTOOL_AMPI_pop_request() { 
                0292   return 0;
                0293 }
                0294 
                0295 void ADTOOL_AMPI_push_AMPI_WinRequest(AMPI_WinRequest *winRequest) {
                0296 }
                0297 
                0298 void ADTOOL_AMPI_pop_AMPI_WinRequest(AMPI_WinRequest *winRequest) {
                0299 }
                0300 
                0301 void ADTOOL_AMPI_push_comm(MPI_Comm comm) {
                0302 }
                0303 
                0304 MPI_Comm ADTOOL_AMPI_pop_comm() {
                0305   return 0;
                0306 }
                0307 
                0308 /**
                0309  * Register the info that the shadow communicator "dupComm"
                0310  * has been created for the new communicator "comm"
                0311  */
                0312 void ADTOOL_AMPI_addShadowComm(MPI_Comm comm, MPI_Comm dupComm) {
                0313   struct AMPI_ShadowComm_list *newCell =
                0314     (struct AMPI_ShadowComm_list *)malloc(sizeof(struct AMPI_ShadowComm_list)) ;
                0315   newCell->next_p = ADTOOL_AMPI_SHADOWCOMM_LIST ;
                0316   newCell->comm = comm;
                0317   newCell->shadowComm = dupComm;
                0318   ADTOOL_AMPI_SHADOWCOMM_LIST = newCell;
                0319 }
                0320 
                0321 /**
                0322  * Get the shadow communicator used to separate the communication graph of
                0323  * (tangent-)diff variables from the communication graph of original variables
                0324  */
                0325 MPI_Comm ADTOOL_AMPI_getShadowComm(MPI_Comm comm) {
                0326   struct AMPI_ShadowComm_list * inShadowCommList = ADTOOL_AMPI_SHADOWCOMM_LIST ;
                0327   while (inShadowCommList!=NULL && inShadowCommList->comm!=comm) {
                0328     inShadowCommList = inShadowCommList->next_p ;
                0329   }
                0330   if (inShadowCommList) {
                0331     return inShadowCommList->shadowComm ;
                0332   } else {
                0333     /* Nothing found about "comm": this is wrong!! fallback return comm */
                0334     return comm ;
                0335   }
                0336 }
                0337 
                0338 /**
                0339  * Removes the info about the shadow communicator associated to "comm".
                0340  */
                0341 void ADTOOL_AMPI_delShadowComm(MPI_Comm comm) {
                0342   struct AMPI_ShadowComm_list ** toinShadowCommList = &ADTOOL_AMPI_SHADOWCOMM_LIST ;
                0343   while (*toinShadowCommList!=NULL && (*toinShadowCommList)->comm!=comm) {
                0344     toinShadowCommList = &((*toinShadowCommList)->next_p) ;
                0345   }
                0346   if (*toinShadowCommList!=NULL) {
                0347     struct AMPI_ShadowComm_list *cell = *toinShadowCommList ;
                0348     toinShadowCommList = &(cell->next_p) ;
                0349     free(cell);
                0350   }
                0351 }
                0352 
                0353 /** Returns the non-diff part of a communication buffer
                0354  * passed to AMPI send or recv. For Tapenade, this is
                0355  * the communication buffer itself (association by name) */
                0356 void* ADTOOL_AMPI_rawData(void* activeData, int *size) { 
                0357   return activeData ;
                0358 }
                0359 
                0360 /**
                0361  * see \ref ADTOOL_AMPI_rawData
                0362  */
                0363 void* ADTOOL_AMPI_rawDataV(void* activeData, int commSize,  int *counts, int* displs) {
                0364   return activeData;
                0365 }
                0366 
                0367 /**
                0368  * returns contiguous data from indata
                0369  */
                0370 void * ADTOOL_AMPI_packDType(void* indata, void* outdata, int count, int idx) {
                0371   return indata;
                0372 }
                0373 
                0374 /**
                0375  * unpacks contiguous data back into structure
                0376  */
                0377 void * ADTOOL_AMPI_unpackDType(void* indata, void* outdata, int count, int idx) {
                0378   return indata;
                0379 }
                0380 
                0381 /** Returns the diff part of the adjoint of a communication buffer
                0382  * passed to AMPI send or recv. For Tapenade, this is
                0383  * the adjoint communication buffer itself (association by name) */
                0384 void* ADTOOL_AMPI_rawAdjointData(void* activeData) { 
                0385   return activeData ;
                0386 }
                0387 
                0388 /** Remembers the association from a request <tt>ampiRequest</tt> to its
                0389  * associated non-diff buffer <tt>buf</tt> */
                0390 void ADTOOL_AMPI_mapBufForAdjoint(struct AMPI_Request_S  *ampiRequest,
                0391                   void* buf) { 
                0392   ampiRequest->buf = buf ;
                0393   ampiRequest->adjointBuf = NULL ;
                0394 }
                0395 
                0396 /** Adds into the request-to-buffer association list the associated
                0397  * adjoint buffer <tt>adjointBuf</tt> of non-diff buffer <tt>buf</tt>
                0398  * This should be done upon turn from FW sweep to BW sweep. */
                0399 void ADTOOL_AMPI_Turn(void* buf, void* adjointBuf) {
                0400   struct AMPI_Request_stack* inStack = requestStackTop ;
                0401   while (inStack!=NULL) {
                0402     if (inStack->buf==buf) {
                0403       inStack->adjointBuf = adjointBuf ;
                0404     }
                0405     inStack = inStack->next_p ;
                0406   }
                0407 }
                0408 
                0409 /*[llh] not used ? redundant with mapBufForAdjoint? */
                0410 void ADTOOL_AMPI_setBufForAdjoint(struct AMPI_Request_S  *ampiRequest,
                0411                   void* buf) { 
                0412   /* an overloading tool would not do this but rather leave the buffer as traced 
                0413      because the memory mapping happens already at FW time */
                0414   ampiRequest->buf=buf;
                0415 }
                0416 
                0417 void ADTOOL_AMPI_getAdjointCount(int *count,
                0418                  MPI_Datatype datatype) {
                0419   /* for now we keep the count as is but for example in vector mode one would have to multiply by vector length */
                0420 }
                0421 
                0422 void ADTOOL_AMPI_setAdjointCount(struct AMPI_Request_S  *ampiRequest) { 
                0423   ampiRequest->adjointCount=ampiRequest->count;
                0424   ADTOOL_AMPI_getAdjointCount(&(ampiRequest->adjointCount),ampiRequest->datatype);
                0425 }
                0426 
                0427 void ADTOOL_AMPI_setAdjointCountAndTempBuf(struct AMPI_Request_S *ampiRequest) { 
                0428   ADTOOL_AMPI_setAdjointCount(ampiRequest);
                0429   ampiRequest->adjointTempBuf =
                0430     ADTOOL_AMPI_allocateTempBuf(ampiRequest->adjointCount,
                0431                                 ampiRequest->datatype,
                0432                                 ampiRequest->comm) ;
                0433   assert(ampiRequest->adjointTempBuf);
                0434 }
                0435 
                0436 void* ADTOOL_AMPI_allocateTempBuf(int adjointCount, MPI_Datatype datatype, MPI_Comm comm) {
                0437   size_t s=0;
                0438   int dt_idx = derivedTypeIdx(datatype);
                0439   if (datatype==MPI_DOUBLE || datatype==MPI_DOUBLE_PRECISION)
                0440     s=sizeof(double);
                0441   else if (datatype==MPI_FLOAT || datatype==MPI_REAL)
                0442     s=sizeof(float);
                0443   else if (isDerivedType(dt_idx))
                0444     s = getDTypeData()->p_extents[dt_idx];
                0445   else
                0446     MPI_Abort(comm, MPI_ERR_TYPE);
                0447   return (void*)malloc(adjointCount*s);
                0448 }
                0449 
                0450 void ADTOOL_AMPI_releaseAdjointTempBuf(void *tempBuf) { 
                0451   free(tempBuf) ;
                0452 }
                0453 
                0454 void* ADTOOL_AMPI_allocateTempActiveBuf(int count,
                0455                     MPI_Datatype datatype,
                0456                     MPI_Comm comm) {
                0457 /*   void* ptr = NULL; */
                0458 /*   if (datatype==MPI_DOUBLE) */
                0459 /*     ptr = malloc(count*sizeof(MPI_DOUBLE)); */
                0460 /*   else if (datatype==MPI_FLOAT) */
                0461 /*     ptr = malloc(count*sizeof(MPI_FLOAT)); */
                0462   MPI_Aint lb, extent ;
                0463   int rc = MPI_Type_get_extent(datatype, &lb, &extent) ;
                0464   assert(rc==MPI_SUCCESS);
                0465   void* ptr = NULL ;
                0466   int size = ((char*)extent)-((char*)lb) ;
                0467   ptr = malloc(count*size) ;
                0468   assert(ptr);
                0469   return ptr;
                0470 }
                0471 
                0472 void ADTOOL_AMPI_releaseTempActiveBuf(void *buf,
                0473                       int count,
                0474                       MPI_Datatype datatype) {
                0475   free(buf);
                0476 }
                0477 
                0478 void *ADTOOL_AMPI_copyActiveBuf(void* source,
                0479                                 void* target,
                0480                                 int count,
                0481                                 MPI_Datatype datatype,
                0482                                 MPI_Comm comm) {
                0483   MPI_Aint lb,extent ;
                0484   int rc = MPI_Type_get_extent(datatype, &lb, &extent) ;
                0485   assert(rc==MPI_SUCCESS);
                0486   int size = extent - lb ;
                0487   memcpy(target, source, count*size) ;
                0488   return source;
                0489 }
                0490 
                0491 /** This is the tangent of assignment target=source*target.
                0492  * If targetd is NULL, does only the original assignment */
                0493 void ADTOOL_AMPI_tangentMultiply(int count, MPI_Datatype datatype, MPI_Comm comm,
                0494                                  void *source, void *tangentSource,
                0495                                  void* target, void* tangentTarget) {
                0496   int i ;
                0497   if (datatype==MPI_DOUBLE || datatype==MPI_DOUBLE_PRECISION) {
                0498     double* tgt = (double*)target ;
                0499     double* tgtd = (double*)tangentTarget ;
                0500     double* src = (double*)source ;
                0501     double* srcd = (double*)tangentSource ;
                0502     for (i=0 ; i<count ; ++i) {
                0503       if (tgtd) tgtd[i] = tgtd[i]*src[i] + tgt[i]*srcd[i] ;
                0504       tgt[i] = tgt[i]*src[i] ;
                0505     }
                0506   } else if (datatype==MPI_FLOAT || datatype==MPI_REAL) {
                0507     float* tgt = (float*)target ;
                0508     float* tgtd = (float*)tangentTarget ;
                0509     float* src = (float*)source ;
                0510     float* srcd = (float*)tangentSource ;
                0511     for (i=0 ; i<count ; ++i) {
                0512       if (tgtd) tgtd[i] = tgtd[i]*src[i] + tgt[i]*srcd[i] ;
                0513       tgt[i] = tgt[i]*src[i] ;
                0514     }
                0515   } else
                0516     MPI_Abort(comm, MPI_ERR_TYPE);
                0517 }
                0518 
                0519 /** This is the tangent of assignment target=MIN(source,target).
                0520  * If targetd is NULL, does only the original assignment */
                0521 void ADTOOL_AMPI_tangentMin(int count, MPI_Datatype datatype, MPI_Comm comm,
                0522                             void *source, void *tangentSource,
                0523                             void* target, void* tangentTarget) {
                0524   int i ;
                0525   if (datatype==MPI_DOUBLE || datatype==MPI_DOUBLE_PRECISION) {
                0526     double* tgt = (double*)target ;
                0527     double* tgtd = (double*)tangentTarget ;
                0528     double* src = (double*)source ;
                0529     double* srcd = (double*)tangentSource ;
                0530     for (i=0 ; i<count ; ++i) {
                0531       if (tgt[i] > src[i]) {
                0532         if (tgtd) tgtd[i] = srcd[i] ;
                0533         tgt[i] = src[i] ;
                0534       }
                0535     }
                0536   } else if (datatype==MPI_FLOAT || datatype==MPI_REAL) {
                0537     float* tgt = (float*)target ;
                0538     float* tgtd = (float*)tangentTarget ;
                0539     float* src = (float*)source ;
                0540     float* srcd = (float*)tangentSource ;
                0541     for (i=0 ; i<count ; ++i) {
                0542       if (tgt[i] > src[i]) {
                0543         if (tgtd) tgtd[i] = srcd[i] ;
                0544         tgt[i] = src[i] ;
                0545       }
                0546     }
                0547   } else
                0548     MPI_Abort(comm, MPI_ERR_TYPE);
                0549 }
                0550 
                0551 /** This is the tangent of assignment target=MAX(source,target).
                0552  * If targetd is NULL, does only the original assignment */
                0553 void ADTOOL_AMPI_tangentMax(int count, MPI_Datatype datatype, MPI_Comm comm,
                0554                             void *source, void *tangentSource,
                0555                             void* target, void* tangentTarget) {
                0556   int i ;
                0557   if (datatype==MPI_DOUBLE || datatype==MPI_DOUBLE_PRECISION) {
                0558     double* tgt = (double*)target ;
                0559     double* tgtd = (double*)tangentTarget ;
                0560     double* src = (double*)source ;
                0561     double* srcd = (double*)tangentSource ;
                0562     for (i=0 ; i<count ; ++i) {
                0563       if (tgt[i] < src[i]) {
                0564         if (tgtd) tgtd[i] = srcd[i] ;
                0565         tgt[i] = src[i] ;
                0566       }
                0567     }
                0568   } else if (datatype==MPI_FLOAT || datatype==MPI_REAL) {
                0569     float* tgt = (float*)target ;
                0570     float* tgtd = (float*)tangentTarget ;
                0571     float* src = (float*)source ;
                0572     float* srcd = (float*)tangentSource ;
                0573     for (i=0 ; i<count ; ++i) {
                0574       if (tgt[i] < src[i]) {
                0575         if (tgtd) tgtd[i] = srcd[i] ;
                0576         tgt[i] = src[i] ;
                0577       }
                0578     }
                0579   } else
                0580     MPI_Abort(comm, MPI_ERR_TYPE);
                0581 }
                0582 
                0583 /**
                0584  * This is the adjoint of assignment target=source*target
                0585  */
                0586 void ADTOOL_AMPI_adjointMultiply(int count, MPI_Datatype datatype, MPI_Comm comm,
                0587                                  void *source, void *adjointSource,
                0588                                  void* target, void* adjointTarget) {
                0589   int i ;
                0590   if (datatype==MPI_DOUBLE || datatype==MPI_DOUBLE_PRECISION) {
                0591     double* tgt = (double*)target ;
                0592     double* tgtb = (double*)adjointTarget ;
                0593     double* src = (double*)source ;
                0594     double* srcb = (double*)adjointSource ;
                0595     for (i=0 ; i<count ; ++i) {
                0596       srcb[i] += tgt[i]*tgtb[i] ;
                0597       tgtb[i] *= src[i] ;
                0598     }
                0599   } else if (datatype==MPI_FLOAT || datatype==MPI_REAL) {
                0600     float* tgt = (float*)target ;
                0601     float* tgtb = (float*)adjointTarget ;
                0602     float* src = (float*)source ;
                0603     float* srcb = (float*)adjointSource ;
                0604     for (i=0 ; i<count ; ++i) {
                0605       srcb[i] += tgt[i]*tgtb[i] ;
                0606       tgtb[i] *= src[i] ;
                0607     }
                0608   } else
                0609     MPI_Abort(comm, MPI_ERR_TYPE);
                0610 }
                0611 
                0612 /**
                0613  * This is the adjoint of assignment target=MIN(source,target)
                0614  */
                0615 void ADTOOL_AMPI_adjointMin(int count, MPI_Datatype datatype, MPI_Comm comm,
                0616                                  void *source, void *adjointSource,
                0617                                  void* target, void* adjointTarget) {
                0618   int i ;
                0619   if (datatype==MPI_DOUBLE || datatype==MPI_DOUBLE_PRECISION) {
                0620     double* tgt = (double*)target ;
                0621     double* tgtb = (double*)adjointTarget ;
                0622     double* src = (double*)source ;
                0623     double* srcb = (double*)adjointSource ;
                0624     for (i=0 ; i<count ; ++i) {
                0625       if (src[i]<tgt[i]) {
                0626         srcb[i] += tgtb[i] ;
                0627         tgtb[i] = 0.0 ;
                0628       }
                0629     }
                0630   } else if (datatype==MPI_FLOAT || datatype==MPI_REAL) {
                0631     float* tgt = (float*)target ;
                0632     float* tgtb = (float*)adjointTarget ;
                0633     float* src = (float*)source ;
                0634     float* srcb = (float*)adjointSource ;
                0635     for (i=0 ; i<count ; ++i) {
                0636       if (src[i]<tgt[i]) {
                0637         srcb[i] += tgtb[i] ;
                0638         tgtb[i] = 0.0 ;
                0639       }
                0640     }
                0641   } else
                0642     MPI_Abort(comm, MPI_ERR_TYPE);
                0643 }
                0644 
                0645 /**
                0646  * This is the adjoint of assignment target=MAX(source,target)
                0647  */
                0648 void ADTOOL_AMPI_adjointMax(int count, MPI_Datatype datatype, MPI_Comm comm,
                0649                                  void *source, void *adjointSource,
                0650                                  void* target, void* adjointTarget) {
                0651   int i ;
                0652   if (datatype==MPI_DOUBLE || datatype==MPI_DOUBLE_PRECISION) {
                0653     double* tgt = (double*)target ;
                0654     double* tgtb = (double*)adjointTarget ;
                0655     double* src = (double*)source ;
                0656     double* srcb = (double*)adjointSource ;
                0657     for (i=0 ; i<count ; ++i) {
                0658       if (src[i]>tgt[i]) {
                0659         srcb[i] += tgtb[i] ;
                0660         tgtb[i] = 0.0 ;
                0661       }
                0662     }
                0663   } else if (datatype==MPI_FLOAT || datatype==MPI_REAL) {
                0664     float* tgt = (float*)target ;
                0665     float* tgtb = (float*)adjointTarget ;
                0666     float* src = (float*)source ;
                0667     float* srcb = (float*)adjointSource ;
                0668     for (i=0 ; i<count ; ++i) {
                0669       if (src[i]>tgt[i]) {
                0670         srcb[i] += tgtb[i] ;
                0671         tgtb[i] = 0.0 ;
                0672       }
                0673     }
                0674   } else
                0675     MPI_Abort(comm, MPI_ERR_TYPE);
                0676 }
                0677 
                0678 /**
                0679  * Multiply the given buffer target, which holds an adjoint, with the given source value
                0680  */
                0681 void ADTOOL_AMPI_multiplyAdjoint(int adjointCount, MPI_Datatype datatype, MPI_Comm comm, void* target, void *source, void *idx) { 
                0682   if (datatype==MPI_DOUBLE || datatype==MPI_DOUBLE_PRECISION) {
                0683     double *vb = (double *)target ;
                0684     double *nb = (double *)source ;
                0685     int i ;
                0686     for (i=0 ; i<adjointCount ; ++i) {
                0687       *vb = *vb * (*nb) ;
                0688       ++vb ;
                0689       ++nb ;
                0690     }
                0691   } else if (datatype==MPI_FLOAT || datatype==MPI_REAL) {
                0692     float *vb = (float *)target ;
                0693     float *nb = (float *)source ;
                0694     int i ;
                0695     for (i=0 ; i<adjointCount ; ++i) {
                0696       *vb = *vb * (*nb) ;
                0697       ++vb ;
                0698       ++nb ;
                0699     }
                0700   } else
                0701     MPI_Abort(comm, MPI_ERR_TYPE);
                0702 }
                0703 
                0704 /**
                0705  * Divide the given buffer target, which holds an adjoint, with the given source value
                0706  */
                0707 void ADTOOL_AMPI_divideAdjoint(int adjointCount, MPI_Datatype datatype, MPI_Comm comm, void* target, void *source, void *idx) { 
                0708   if (datatype==MPI_DOUBLE || datatype==MPI_DOUBLE_PRECISION) {
                0709     double *vb = (double *)target ;
                0710     double *nb = (double *)source ;
                0711     int i ;
                0712     for (i=0 ; i<adjointCount ; ++i) {
                0713       *vb = *vb / *nb ;
                0714       ++vb ;
                0715       ++nb ;
                0716     }
                0717   } else if (datatype==MPI_FLOAT || datatype==MPI_REAL) {
                0718     float *vb = (float *)target ;
                0719     float *nb = (float *)source ;
                0720     int i ;
                0721     for (i=0 ; i<adjointCount ; ++i) {
                0722       *vb = *vb / *nb ;
                0723       ++vb ;
                0724       ++nb ;
                0725     }
                0726   } else
                0727     MPI_Abort(comm, MPI_ERR_TYPE);
                0728 }
                0729 
                0730 /**
                0731  * Check equality of the given buffers source1 and source2, which hold adjoints,
                0732  * and return the result in the given target buffer.
                0733  */
                0734 void ADTOOL_AMPI_equalAdjoints(int adjointCount, MPI_Datatype datatype, MPI_Comm comm, void* target, void *source1, void *source2, void *idx) { 
                0735   if (datatype==MPI_DOUBLE || datatype==MPI_DOUBLE_PRECISION) {
                0736     double *vb = (double *)target ;
                0737     double *nb = (double *)source1 ;
                0738     double *fb = (double *)source2 ;
                0739     int i ;
                0740     for (i=0 ; i<adjointCount ; ++i) {
                0741       *vb = *nb == *fb ;
                0742       ++vb ;
                0743       ++nb ;
                0744       ++fb ;
                0745     }
                0746   } else if (datatype==MPI_FLOAT || datatype==MPI_REAL) {
                0747     float *vb = (float *)target ;
                0748     float *nb = (float *)source1 ;
                0749     float *fb = (float *)source2 ;
                0750     int i ;
                0751     for (i=0 ; i<adjointCount ; ++i) {
                0752       *vb = *nb == *fb ;
                0753       ++vb ;
                0754       ++nb ;
                0755       ++fb ;
                0756     }
                0757   } else
                0758     MPI_Abort(comm, MPI_ERR_TYPE);
                0759 }
                0760 
                0761 /**
                0762  * Increment the given buffer "target", which holds an adjoint variable,
                0763  * with the given additional adjoint value found in "source".
                0764  */
                0765 void ADTOOL_AMPI_incrementAdjoint(int adjointCount, MPI_Datatype datatype, MPI_Comm comm, void* target, void *source, void *idx) { 
                0766   int dt_idx = derivedTypeIdx(datatype);
                0767   if (isUserDefinedOp(dt_idx)) {
                0768     derivedTypeData* dat = getDTypeData();
                0769     MPI_Aint lb, extent ;
                0770     MPI_Type_get_extent(datatype,&lb,&extent);
                0771     MPI_Aint*   fieldOffsets = dat->arrays_of_displacements[dt_idx] ;
                0772     int*   fieldBlocklengths = dat->arrays_of_blocklengths[dt_idx] ;
                0773     MPI_Datatype* fieldTypes = dat->arrays_of_types[dt_idx] ;
                0774     int nbfields = dat->counts[dt_idx] ;
                0775     int i,j ;
                0776     for (i=0 ; i<adjointCount ; ++i) {
                0777       for (j=0 ; j<nbfields ; ++j) {
                0778         ADTOOL_AMPI_incrementAdjoint(fieldBlocklengths[j], fieldTypes[j], comm,
                0779                                      target+fieldOffsets[j], source+fieldOffsets[j], idx) ;
                0780       }
                0781       target += extent ;
                0782       source += extent ;
                0783     }
                0784   } else if (datatype==MPI_DOUBLE || datatype==MPI_DOUBLE_PRECISION) {
                0785     double *vb = (double *)target ;
                0786     double *nb = (double *)source ;
                0787     int i ;
                0788     for (i=0 ; i<adjointCount ; ++i) {
                0789       *vb = *vb + *nb ;
                0790       ++vb ;
                0791       ++nb ;
                0792     }
                0793   } else if (datatype==MPI_FLOAT || datatype==MPI_REAL) {
                0794     float *vb = (float *)target ;
                0795     float *nb = (float *)source ;
                0796     int i ;
                0797     for (i=0 ; i<adjointCount ; ++i) {
                0798       *vb = *vb + *nb ;
                0799       ++vb ;
                0800       ++nb ;
                0801     }
                0802   } else
                0803     MPI_Abort(comm, MPI_ERR_TYPE);
                0804 }
                0805 
                0806 /**
                0807  * Reset to zero the given buffer "target", which holds an adjoint variable.
                0808  */
                0809 void ADTOOL_AMPI_nullifyAdjoint(int adjointCount, MPI_Datatype datatype, MPI_Comm comm,
                0810                                 void* target) {
                0811   int dt_idx = derivedTypeIdx(datatype);
                0812   if (isUserDefinedOp(dt_idx)) {
                0813     derivedTypeData* dat = getDTypeData();
                0814     MPI_Aint lb, extent ;
                0815     MPI_Type_get_extent(datatype,&lb,&extent);
                0816     MPI_Aint*   fieldOffsets = dat->arrays_of_displacements[dt_idx] ;
                0817     int*   fieldBlocklengths = dat->arrays_of_blocklengths[dt_idx] ;
                0818     MPI_Datatype* fieldTypes = dat->arrays_of_types[dt_idx] ;
                0819     int nbfields = dat->counts[dt_idx] ;
                0820     int i,j ;
                0821     for (i=0 ; i<adjointCount ; ++i) {
                0822       for (j=0 ; j<nbfields ; ++j) {
                0823         ADTOOL_AMPI_nullifyAdjoint(fieldBlocklengths[j], fieldTypes[j], comm,
                0824                                    target+fieldOffsets[j]) ;
                0825       }
                0826       target += extent ;
                0827     }
                0828   } else if (datatype==MPI_DOUBLE || datatype==MPI_DOUBLE_PRECISION) {
                0829     double *vb = (double *)target ;
                0830     int i ;
                0831     for (i=0 ; i<adjointCount ; ++i) {
                0832       *vb = 0.0 ;
                0833       ++vb ;
                0834     }
                0835   } else if (datatype==MPI_FLOAT || datatype==MPI_REAL) {
                0836     float *vb = (float *)target ;
                0837     int i ;
                0838     for (i=0 ; i<adjointCount ; ++i) {
                0839       *vb = 0.0 ;
                0840       ++vb ;
                0841     }
                0842   } else
                0843     MPI_Abort(comm, MPI_ERR_TYPE);
                0844 }
                0845 
                0846 extern void pushNArray(void *x, unsigned int nbChars) ;
                0847 extern void popNArray(void *x, unsigned int nbChars) ;
                0848 
                0849 /**
                0850  * Push the contents of buffer somewhere
                0851  */
                0852 void ADTOOL_AMPI_pushBuffer(int count, MPI_Datatype datatype, MPI_Comm comm,
                0853                             void* buffer) {
                0854     MPI_Aint lb, extent ;
                0855     MPI_Type_get_extent(datatype,&lb,&extent);
                0856     int length = count*(int)extent ;
                0857     pushNArray((char*)buffer, length) ;
                0858 }
                0859 
                0860 /**
                0861  * Pop the contents of buffer from somewhere
                0862  */
                0863 void ADTOOL_AMPI_popBuffer(int count, MPI_Datatype datatype, MPI_Comm comm,
                0864                            void* buffer) {
                0865     MPI_Aint lb, extent ;
                0866     MPI_Type_get_extent(datatype,&lb,&extent);
                0867     int length = count*(int)extent ;
                0868     popNArray((char*)buffer, length) ;
                0869 }
                0870 
                0871 void ADTOOL_AMPI_writeData(void *buf,int *count) {}
                0872 
                0873 void ADTOOL_AMPI_writeDataV(void* activeData, int *counts, int* displs) {}
                0874 
                0875 void ADTOOL_AMPI_setupTypes() {
                0876 #ifdef AMPI_FORTRANCOMPATIBLE
                0877   MPI_Fint adouble;
                0878   MPI_Fint areal;
                0879 #endif
                0880   /* Change AMPI_ADOUBLE to something else? Need AMPI_ADOUBLE!=MPI_DOUBLE for derived types. */
                0881   AMPI_ADOUBLE=MPI_DOUBLE;
                0882   AMPI_AFLOAT=MPI_FLOAT;
                0883 #ifdef AMPI_FORTRANCOMPATIBLE
                0884   adtool_ampi_fortransetuptypes_(&adouble, &areal);
                0885   AMPI_ADOUBLE_PRECISION=MPI_Type_f2c(adouble);
                0886   AMPI_AREAL=MPI_Type_f2c(areal);
                0887 #endif
                0888 }
                0889 
                0890 void ADTOOL_AMPI_cleanupTypes() {
                0891 #ifdef AMPI_FORTRANCOMPATIBLE
                0892   MPI_Fint adouble=MPI_Type_c2f(AMPI_ADOUBLE_PRECISION);
                0893   MPI_Fint areal=MPI_Type_c2f(AMPI_AREAL);
                0894   adtool_ampi_fortrancleanuptypes_(&adouble, &areal);
                0895 #endif
                0896 }
                0897 
                0898 MPI_Datatype ADTOOL_AMPI_FW_rawType(MPI_Datatype datatype) {
                0899   int dt_idx = derivedTypeIdx(datatype);
                0900   if (datatype==AMPI_ADOUBLE) return MPI_DOUBLE;
                0901   else if (datatype==AMPI_AFLOAT) return MPI_FLOAT;
                0902   else if (isDerivedType(dt_idx)) return getDTypeData()->packed_types[dt_idx];
                0903   else return datatype;
                0904 }
                0905 
                0906 MPI_Datatype ADTOOL_AMPI_BW_rawType(MPI_Datatype datatype) {
                0907   int dt_idx = derivedTypeIdx(datatype);
                0908   if (datatype==AMPI_ADOUBLE) return MPI_DOUBLE;
                0909   else if (datatype==AMPI_AFLOAT) return MPI_FLOAT;
                0910   else if (isDerivedType(dt_idx)) return MPI_DOUBLE;
                0911   else return datatype;
                0912 }
                0913 
                0914 AMPI_Activity ADTOOL_AMPI_isActiveType(MPI_Datatype datatype) {
                0915   if (datatype==AMPI_ADOUBLE
                0916       ||
                0917       datatype==AMPI_AFLOAT
                0918 #ifdef AMPI_FORTRANCOMPATIBLE
                0919       ||
                0920       datatype==AMPI_ADOUBLE_PRECISION
                0921       ||
                0922       datatype==AMPI_AREAL
                0923 #endif
                0924       ) return AMPI_ACTIVE;
                0925   return AMPI_PASSIVE;
                0926 }
                0927 
                0928 void *ADTOOL_AMPI_createWinMap(void *active_buf, MPI_Aint size){
                0929   return NULL;
                0930 }
                0931 
                0932 void ADTOOL_AMPI_writeWinData(void *map, void *buf, MPI_Aint size){
                0933 }
                0934 
                0935 MPI_Aint ADTOOL_AMPI_getWinSize(MPI_Aint size) {
                0936    return 0;
                0937 }